Bayesian inference

Inference for age structured compartment models using the diffusion approximation (via the van Kampen expansion). See this paper for more details on the method.

There are two ways to do inference: manifold method (sec 3.3 in the report) and tangent space method (sec 3.4 in the report). In various degrees of less robust but fast to more robust but slow:

  • tangent space method.
  • manifold method with few internal steps and fast integration method (det_method = RK2, lyapunov_method = euler).
  • manifold method with large number of internel steps and robust integration method (solve_ivp from scipy library).
Methods for full data  
infer Infers epidemiological and control parameters given all information.
infer_mcmc Explore the posterior distribution given all information.
infer_nested_sampling Compute the model evidence (and generate posterior samples) given all information.
obtain_minus_log_p Computes -log(p) of a fully observed trajectory.
compute_hessian Computes the Hessian of -log(p).
nested_sampling_inference Compute the log-evidence and weighted samples.
Methods for partial data  
latent_infer Infers epidemiological and control parameters and initial conditions.
latent_infer_mcmc Explore the posterior distribution.
latent_infer_nested_sampling Compute the model evidence (and generate posterior samples).
minus_logp_red Computes -log(p) of a partially observed trajectory.
compute_hessian_latent Computes the Hessian of -log(p).
nested_sampling_latent_inference Compute the log-evidence and weighted samples.
Sensitivity analysis  
FIM Computes the Fisher Information Matrix of the stochastic model.
FIM_det Computes the Fisher Information Matrix of the deterministic model.
sensitivity Computes the normalized sensitivity measure
Helper function  
integrate A wrapper around ‘simulate’ in pyross.deterministic.
set_params Sets parameters.
set_det_method Sets the integration method of the deterministic equation
set_lyapunov_method Sets the integration method of the Lyapunov equation
set_det_model Sets the internal deterministic model
set_contact_matrix Sets the contact matrix
fill_params_dict Fills and returns a parameter dictionary
get_mean_inits Constructs full initial conditions from the prior dict

The functions are documented under the parent class SIR_type.

SIR_type

class pyross.inference.SIR_type

Parent class for inference for all SIR-type classes listed below

All subclasses use the same functions to perform inference, which are documented below.

FIM()

Computes the Fisher Information Matrix (FIM) for the MAP estimates of a stochastic SIR type model.

Parameters:
  • x (2d numpy.array) – Observed trajectory (number of data points x (age groups * model classes))
  • Tf (float) – Total time of the trajectory
  • infer_result (dict) – Dictionary returned by latent_infer
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • eps (float or numpy.array, optional) –

    Step size for numerical differentiation of the process mean and its full covariance matrix with respect to the parameters. Must be either a scalar, or an array of length len(infer_result[‘flat_params’]). If not specified,

    eps = 100*infer_result['flat_params']
          *numpy.divide(numpy.spacing(infer_result['log_likelihood']),
          infer_result['log_likelihood'])**(0.25)
    

    is used. It is recommended to use a step-size greater or equal to eps. Decreasing the step size too small can result in round-off error.

  • inter_steps (int, optional) – Intermediate steps for interpolation between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration. We have found that forward Euler is generally slower, but more stable for derivatives with respect to parameters than the variable step size integrators used elsewhere in pyross. Default is 100.
Returns:

FIM – The Fisher Information Matrix

Return type:

2d numpy.array

FIM_det()

Computes the Fisher Information Matrix (FIM) for the MAP estimates of a deterministic (ODE based, including a constant measurement error) SIR type model.

Parameters:
  • x (2d numpy.array) – Observed trajectory (number of data points x (age groups * model classes))
  • Tf (float) – Total time of the trajectory
  • infer_result (dict) – Dictionary returned by latent_infer
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • eps (float or numpy.array, optional) –

    Step size for numerical differentiation of the process mean and its full covariance matrix with respect to the parameters. Must be either a scalar, or an array of length len(infer_result[‘flat_params’]). If not specified,

    eps = 100*infer_result['flat_params']
          *numpy.divide(numpy.spacing(infer_result['log_likelihood']),
          infer_result['log_likelihood'])**(0.25)
    

    is used. It is recommended to use a step-size greater or equal to eps. Decreasing the step size too small can result in round-off error.

  • measurement_error (float, optional) – Standard deviation of measurements (uniform and independent Gaussian measurement error assumed). Default is 1e-2.
  • inter_steps (int, optional) – Intermediate steps for interpolation between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration. We have found that forward Euler is generally slower, but more stable for derivatives with respect to parameters than the variable step size integrators used elsewhere in pyross. Default is 100.
Returns:

FIM – The Fisher Information Matrix

Return type:

2d numpy.array

evidence_laplace()

Compute the evidence using a Laplace approximation at the MAP estimate.

Parameters:
  • x (2d numpy.array) – Observed trajectory (number of data points x (age groups * model classes))
  • Tf (float) – Total time of the trajectory
  • infer_result (dict) – Dictionary returned by latent_infer
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • eps (float or numpy.array, optional) –

    Step size for finite differences computation of the hessian with respect to the parameters. Must be either a scalar, or an array of length len(infer_result[‘flat_params’]). If not specified,

    eps = 100*infer_result['flat_params']
          *numpy.divide(numpy.spacing(infer_result['log_likelihood']),
          infer_result['log_likelihood'])**(0.25)
    

    is used. For fd_method=”central” it is recommended to use a step-size greater or equal to eps. Decreasing the step size too small can result in round-off error.

  • fd_method (str, optional) – The type of finite-difference scheme used to compute the hessian, supports “forward” and “central”. Default is “central”.
  • inter_steps (int, optional) – Only used if tangent=False. Intermediate steps for interpolation between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration. We have found that forward Euler is generally slower, but more stable for derivatives with respect to parameters than the variable step size integrators used elsewhere in pyross. Default is 10.
Returns:

log_evidence – The log-evidence computed via Laplace approximation at the MAP estimate.

Return type:

float

fill_params_dict()

Returns a full dictionary for epidemiological parameters with some changed values

Parameters:
  • keys (list of String) – A list of names of parameters to be changed.
  • params (numpy.array of list) – An array of the same size as keys for the updated value.
  • return_additional_params (boolean, optional (default = False)) – Handling of parameters that are not model parameters (e.g. control parameters). False: raise exception, True: return second dictionary with other parameters
Returns:

full_parameters – A dictionary of epidemiological parameters. For parameter names specified in keys, set the values to be the ones in params; for the others, use the values stored in the class.

Return type:

dict

get_mean_inits()

Construct full initial conditions from the prior dict

Parameters:
  • init_priors (dict) – A dictionary for priors for initial conditions. Same as the init_priors passed to latent_infer. In this function, only takes the mean.
  • obs0 (numpy.array) – Observed initial conditions.
  • fltr0 (numpy.array) – Filter for the observed initial conditions.
Returns:

x0 – Full initial conditions.

Return type:

numpy.array

hessian()

Computes the Hessian matrix for the MAP estimates of an SIR type model.

Parameters:
  • x (2d numpy.array) – Observed trajectory (number of data points x (age groups * model classes))
  • Tf (float) – Total time of the trajectory
  • infer_result (dict) – Dictionary returned by latent_infer
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • eps (float or numpy.array, optional) –

    Step size for finite differences computation of the hessian with respect to the parameters. Must be either a scalar, or an array of length len(infer_result[‘flat_params’]). If not specified,

    eps = 100*infer_result['flat_params']
          *numpy.divide(numpy.spacing(infer_result['log_likelihood']),
          infer_result['log_likelihood'])**(0.25)
    

    is used. For fd_method=”central” it is recommended to use a step-size greater or equal to eps. Decreasing the step size too small can result in round-off error.

  • fd_method (str, optional) – The type of finite-difference scheme used to compute the hessian, supports “forward” and “central”. Default is “central”.
  • inter_steps (int, optional) – Only used if tangent=False. Intermediate steps for interpolation between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration. We have found that forward Euler is generally slower, but sometimes more stable for derivatives with respect to parameters than the variable step size integrators used elsewhere in pyross. Default is 0.
Returns:

hess – The Hessian matrix

Return type:

2d numpy.array

infer()

Compute the maximum a-posteriori (MAP) estimate for all desired parameters, including control parameters, for an SIR type model with fully observed classes. If generator is specified, the lockdown is modelled by scaling the contact matrices for contact at work, school, and other (but not home). This function infers the scaling parameters (can be age dependent) assuming that full data on all classes is available (with latent variables, use latent_infer).

Parameters:
  • x (2d numpy.array) – Observed trajectory (number of data points x (age groups * model classes))
  • Tf (float) – Total time of the trajectory
  • prior_dict (dict) – A dictionary containing priors for parameters (can include both model and intervention parameters). See examples.
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is false.
  • verbose (bool, optional) – Set to True to see intermediate outputs from the optimizer.
  • ftol (double) – Relative tolerance of logp
  • global_max_iter (int, optional) – Number of global optimisations performed.
  • local_max_iter (int, optional) – Number of local optimisation performed.
  • global_atol (float) – The absolute tolerance for global optimisation.
  • enable_global (bool, optional) – Set to True to enable global optimisation.
  • enable_local (bool, optional) – Set to True to enable local optimisation.
  • cma_processes (int, optional) – Number of parallel processes used for global optimisation.
  • cma_population (int, optional) – The number of samples used in each step of the CMA algorithm.
  • cma_random_seed (int (between 0 and 2**32-1)) – Random seed for the optimisation algorithms. By default it is generated from numpy.random.randint.
Returns:

output_dict – Dictionary of MAP estimates, containing the following keys for users:

params_dict: dict

Dictionary for MAP estimates of the model parameters.

control_params_dict: dict

Dictionary for MAP estimates of the control parameters (if requested).

-logp: float

Value of -logp at MAP.

Return type:

dict

Note

This function combines the functionality of infer_parameters and infer_control, which will be deprecated. To infer model parameters only, specify a fixed contactMatrix function. To infer control parameters only, specify a generator and do not specify priors for model parameters.

Examples

An example of prior_dict to set priors for alpha and beta, where alpha is age dependent and we want to infer its scale parameters rather than each component individually. The prior distribution is assumed to be log-normal with the specified mean and standard deviation.

>>> prior_dict = {
        'alpha':{
            'mean': [0.5, 0.2],
            'infer_scale': True,
            'scale_factor_std': 1,
            'scale_factor_bounds': [0.1, 10],
            'prior_fun': 'truncnorm'
        },
        'beta':{
            'mean': 0.02,
            'std': 0.1,
            'bounds': [1e-4, 1],
            'prior_fun': 'lognorm'
        }
    }
infer_control()

Compute the maximum a-posteriori (MAP) estimate of the change of control parameters for a SIR type model in lockdown.

Parameters:infer (see) –
Returns:output_dict – Dictionary of MAP estimates, containing the following keys for users:
map_dict: dict
Dictionary for MAP estimates of the control parameters.
-logp: float
Value of -logp at MAP.
Return type:dict

Note

This function just calls infer with the specified generator, will be deprecated.

infer_mcmc()

Sample the posterior distribution of the epidemiological parameters using ensemble MCMC.

Parameters:
  • x (2d numpy.array) – Observed trajectory (number of data points x (age groups * model classes))
  • Tf (float) – Total time of the trajectory
  • prior_dict (dict) – A dictionary containing priors for parameters (can include both model and intervention parameters). See examples.
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • verbose (bool, optional) – Set to True to see a progress bar for the sample generation. Default is False.
  • sampler (emcee.EnsembleSampler, optional) – Set to instance of the sampler (as returned by this function) to continue running the MCMC chains. Default is None (i.e. run a new chain).
  • nwalkers (int, optional) – The number of chains in the ensemble (should be at least 2*dim). Default is 2*dim.
  • walker_pos (np.array, optional) – The initial position of the walkers. If not specified, it samples random positions from the prior.
  • nsamples (int, optional) – The number of samples per walker. Default is 1000.
  • nprocesses (int, optional) – The number of processes used to compute the likelihood for the walkers, needs pathos. Default is the number of cpu cores if pathos is available, otherwise 1.
Returns:

sampler – This function returns the interal state of the sampler. To look at the chain of the internal flattened parameters, run sampler.get_chain(). Use this to judge whether the chain has sufficiently converged. Either rerun mcmc_inference(…, sampler=sampler) to continue the chain or mcmc_inference_process_result(…) to process the result.

Return type:

emcee.EnsembleSampler

Examples

For the structure of prior_dict, see the documentation of infer. To start sampling the posterior, run for example

>>> sampler = estimator.infer_mcmc(x, Tf, prior_dict, contactMatrix=contactMatrix, verbose=True)

To judge the convergence of this chain, we can look at the trace plot of all the chains (for a moderate number of dimensions dim)

>>> fig, axes = plt.subplots(dim, sharex=True)
>>> samples = sampler.get_chain()
>>> for i in range(dim):
        ax = axes[i]
        ax.plot(samples[:, :, i], "k", alpha=0.3)
        ax.set_xlim(0, len(samples))
>>> axes[-1].set_xlabel("step number");

For more detailed convergence metrics, see the documentation of emcee. To continue running this chain, we can call this function again with the sampler as argument

>>> sampler = estimator.infer_mcmc(x, Tf, prior_dict, contactMatrix=contactMatrix, verbose=True, sampler=sampler)

This procudes 1000 additional samples in each chain. To process the results, call infer_mcmc_process_result.

infer_mcmc_process_result()

Take the sampler generated by pyross.inference.infer_mcmc and produce output dictionaries for further use in the pyross framework. See pyross.inference.infer_mcmc for additional description of parameters.

Parameters:
  • sampler (emcee.EnsembleSampler) – Output of pyross.inference.infer_mcmc.
  • prior_dict (dict) –
  • contactMatrix (callable, optional) –
  • generator (pyross.contactMatrix, optional) –
  • intervention_fun (callable, optional) –
  • flat (bool, optional) – This decides whether to return the samples as for each chain separately (False) or as as a combined list (True). Default is True.
  • discard (int, optional) – The number of initial samples to discard in each chain (to account for burn-in). Default is 0.
  • thin (int, optional) – Thin out the chain by taking only the n-tn element in each chain. Default is 1 (no thinning).
  • **catchall_kwargs (dict) – Catched further provided arguments and ignores them.
Returns:

output_samples – The processed posterior samples.

Return type:

list of dict (if flat=True), or list of list of dict (if flat=False)

infer_nested_sampling()

Compute the log-evidence and weighted samples of the a-posteriori distribution of the parameters of a SIR type model using nested sampling as implemented in the dynesty Python package. This function assumes that full data on all classes is available.

Parameters:
  • x (2d numpy.array) – Observed trajectory (number of data points x (age groups * model classes))
  • Tf (float) – Total time of the trajectory
  • prior_dict (dict) – A dictionary containing priors for parameters (can include both model and intervention parameters). See examples.
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is false.
  • verbose (bool, optional) – Set to True to see intermediate outputs from the nested sampling procedure.
  • nprocesses (int, optional) – The number of processes used for parallel evaluation of the likelihood.
  • queue_size (int, optional) – Size of internal queue of likelihood values, default is nprocesses if multiprocessing is used.
  • maxiter (int, optional) – The maximum number of iterations. Default is no limit.
  • maxcall (int, optional) – The maximum number of calls to the likelihood function. Default no limit.
  • dlogz (float, optional) – The iteration terminates if the estimated contribution of the remaining prior volume to the total evidence falls below this threshold. Default value is 1e-3 * (nlive - 1) + 0.01 if add_live==True, 0.01 otherwise.
  • n_effective (float, optional) – The iteration terminates if the number of effective posterior samples reaches this values. Default is no limit.
  • add_live (bool, optional) – Determines whether to add the remaining set of live points to the set of samples. Default is True.
  • sampler (dynesty.NestedSampler, optional) – Continue running an instance of a nested sampler until the termination criteria are met.
  • **dynesty_args

    Arguments passed through to the construction of the dynesty.NestedSampler constructor. Relevant entries are (this is not comprehensive, for details see the documentation of dynesty):

    nlive: int, optional
    The number of live points. Default is 500.
    bound: {‘none’, ‘single’, ‘multi’, ‘balls’, ‘cubes’}, optional
    Method used to approximately bound the prior using the current set of live points. Default is ‘multi’.
    sample: {‘auto’, ‘unif’, ‘rwalk’, ‘rstagger’, ‘slice’, ‘rslice’, ‘hslice’, callable}, optional
    Method used to sample uniformly within the likelihood constraint, conditioned on the provided bounds.
Returns:

sampler – The state of the sampler after termination of the nested sampling run.

Return type:

dynesty.NestedSampler

infer_nested_sampling_process_result()

Take the sampler generated by pyross.inference.infer_nested_sampling and produce output dictionaries for further use in the pyross framework. See pyross.inference.infer_nested_sampling for description of parameters.

Parameters:
  • sampler (dynesty.NestedSampler) – The output of pyross.inference.infer_nested_sampling.
  • prior_dict (dict) –
  • contactMatrix (callable, optional) –
  • generator (pyross.contactMatrix, optional) –
  • intervention_fun (callable, optional) –
  • **catchall_kwargs (dict) – Catched further provided arguments and ignores them.
Returns:

  • result (dynesty.Result) – The result of the nested sampling iteration. Relevant entries include:

    result.logz: list

    The progression of log-evidence estimates, use result.logz[-1] for the final estimate.

  • output_samples (list) – The processed weighted posterior samples.

infer_parameters()

Infers the MAP estimates for epidemiological parameters

Parameters:infer (see) –
Returns:output – Contains the following keys for users:
map_dict: dict
A dictionary for MAPs. Keys are the names of the parameters and the corresponding values are its MAP estimates.
-logp: float
The value of -logp at MAP.
Return type:dict

Note

This function just calls infer with a fixed contactMatrix function, will be deprecated.

integrate()

An light weight integrate method similar to simulate in pyross.deterministic

Parameters:
  • x0 (np.array) – Initial state of the given model
  • t1 (float) – Initial time of integrator
  • t2 (float) – Final time of integrator
  • steps (int) – Number of time steps for numerical integrator evaluation.
  • max_step (int, optional) – The maximum allowed step size of the integrator.
Returns:

sol – The state of the system evaulated at the time point specified. Only used if det_method is set to ‘solve_ivp’.

Return type:

np.array

latent_FIM()

Computes the Fisher Information Matrix (FIM) of the stochastic model for the initial conditions and all desired parameters, including control parameters, for a SIR type model with partially observed classes. The unobserved classes are treated as latent variables.

Parameters:
  • obs (np.array) – The partially observed trajectory.
  • fltr (2d np.array) – The filter for the observation such that \(F_{ij} x_j (t) = obs_i(t)\)
  • Tf (float) – Total time of the trajectory
  • infer_result (dict) – Dictionary returned by latent_infer
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • eps (float or numpy.array, optional) –

    Step size for numerical differentiation of the process mean and its full covariance matrix with respect to the parameters. Must be either a scalar, or an array of length len(infer_result[‘flat_params’]). If not specified,

    eps = 100*infer_result['flat_params']
          *numpy.divide(numpy.spacing(infer_result['log_likelihood']),
          infer_result['log_likelihood'])**(0.25)
    

    is used. It is recommended to use a step-size greater or equal to eps. Decreasing the step size too small can result in round-off error.

  • inter_steps (int, optional) – Intermediate steps between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration. We have found that forward Euler is generally slower, but more stable for derivatives with respect to parameters than the variable step size integrators used elsewhere in pyross. Default is 100.
Returns:

FIM – The Fisher Information Matrix

Return type:

2d numpy.array

latent_FIM_det()

Computes the Fisher Information Matrix (FIM) of the deterministic model (ODE based, including a constant measurement error) for the initial conditions and all desired parameters, including control parameters, for a SIR type model with partially observed classes. The unobserved classes are treated as latent variables.

Parameters:
  • obs (np.array) – The partially observed trajectory.
  • fltr (2d np.array) – The filter for the observation such that \(F_{ij} x_j (t) = obs_i(t)\)
  • Tf (float) – Total time of the trajectory
  • infer_result (dict) – Dictionary returned by latent_infer
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • eps (float or numpy.array, optional) –

    Step size for numerical differentiation of the process mean and its full covariance matrix with respect to the parameters. Must be either a scalar, or an array of length len(infer_result[‘flat_params’]). If not specified,

    eps = 100*infer_result['flat_params']
          *numpy.divide(numpy.spacing(infer_result['log_likelihood']),
          infer_result['log_likelihood'])**(0.25)
    

    is used. It is recommended to use a step-size greater or equal to eps. Decreasing the step size too small can result in round-off error.

  • measurement_error (float, optional) – Standard deviation of measurements (uniform and independent Gaussian measurement error assumed). Default is 1e-2.
  • inter_steps (int, optional) – Intermediate steps between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration. We have found that forward Euler is generally slower, but more stable for derivatives with respect to parameters than the variable step size integrators used elsewhere in pyross. Default is 100.
Returns:

FIM_det – The Fisher Information Matrix

Return type:

2d numpy.array

latent_evidence_laplace()

Compute the evidence using a Laplace approximation at the MAP estimate for a SIR type model with partially observed classes. The unobserved classes are treated as latent variables.

Parameters:
  • obs (np.array) – The partially observed trajectory.
  • fltr (2d np.array) – The filter for the observation such that \(F_{ij} x_j (t) = obs_i(t)\)
  • Tf (float) – Total time of the trajectory
  • infer_result (dict) – Dictionary returned by latent_infer
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • eps (float or numpy.array, optional) –

    Step size for finite differences computation of the hessian with respect to the parameters. Must be either a scalar, or an array of length len(infer_result[‘flat_params’]). If not specified,

    eps = 100*infer_result['flat_params']
          *numpy.divide(numpy.spacing(infer_result['log_likelihood']),
          infer_result['log_likelihood'])**(0.25)
    

    is used. For fd_method=”central” it is recommended to use a step-size greater or equal to eps. Decreasing the step size too small can result in round-off error.

  • fd_method (str, optional) – The type of finite-difference scheme used to compute the hessian, supports “forward” and “central”. Default is “central”.
  • inter_steps (int, optional) – Intermediate steps between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration. We have found that forward Euler is generally slower, but more stable for derivatives with respect to parameters than the variable step size integrators used elsewhere in pyross. Default is 100.
Returns:

log_evidence – The log-evidence computed via Laplace approximation at the MAP estimate.

Return type:

float

latent_hessian()

Computes the Hessian matrix for the initial conditions and all desired parameters, including control parameters, for a SIR type model with partially observed classes. The unobserved classes are treated as latent variables.

Parameters:
  • obs (np.array) – The partially observed trajectory.
  • fltr (2d np.array) – The filter for the observation such that \(F_{ij} x_j (t) = obs_i(t)\)
  • Tf (float) – Total time of the trajectory
  • infer_result (dict) – Dictionary returned by latent_infer
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • eps (float or numpy.array, optional) –

    Step size for finite differences computation of the hessian with respect to the parameters. Must be either a scalar, or an array of length len(infer_result[‘flat_params’]). If not specified,

    eps = 100*infer_result['flat_params']
          *numpy.divide(numpy.spacing(infer_result['log_likelihood']),
          infer_result['log_likelihood'])**(0.25)
    

    is used. For fd_method=”central” it is recommended to use a step-size greater or equal to eps. Decreasing the step size too small can result in round-off error.

  • fd_method (str, optional) – The type of finite-difference scheme used to compute the hessian, supports “forward” and “central”. Default is “central”.
  • inter_steps (int, optional) – Intermediate steps between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration. We have found that forward Euler is generally slower, but sometimes more stable for derivatives with respect to parameters than the variable step size integrators used elsewhere in pyross. Default is 0.
Returns:

hess – The Hessian matrix

Return type:

2d numpy.array

latent_infer()

Compute the maximum a-posteriori (MAP) estimate for the initial conditions and all desired parameters, including control parameters, for a SIR type model with partially observed classes. The unobserved classes are treated as latent variables.

Parameters:
  • obs (np.array) – The partially observed trajectory.
  • fltr (2d np.array) – The filter for the observation such that \(F_{ij} x_j (t) = obs_i(t)\)
  • Tf (float) – Total time of the trajectory
  • param_priors (dict) – A dictionary that specifies priors for parameters (including control parameters, if desired). See infer for further explanations.
  • init_priors (dict) – A dictionary for priors for initial conditions. See below for examples
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is false.
  • verbose (bool, optional) – Set to True to see intermediate outputs from the optimizer.
  • ftol (double) – Relative tolerance of logp
  • global_max_iter (int, optional) – Number of global optimisations performed.
  • local_max_iter (int, optional) – Number of local optimisation performed.
  • local_initital_step (optional, float or np.array) – Initial step size for the local optimiser. If scalar, relative to the initial guess. Default: Deterined by final state of global optimiser, or, if enable_global=False, 0.01
  • global_atol (float) – The absolute tolerance for global minimisation.
  • enable_global (bool, optional) – Set to True to enable global optimisation.
  • enable_local (bool, optional) – Set to True to enable local optimisation.
  • cma_processes (int, optional) – Number of parallel processes used for global optimisation.
  • cma_population (int, optional) – The number of samples used in each step of the CMA algorithm.
  • cma_random_seed (int (between 0 and 2**32-1)) – Random seed for the optimisation algorithms. By default it is generated from numpy.random.randint.
  • objective (string, optional) – Objective for the minimisation. ‘likelihood’ (default), ‘least_square’ (least squares fit w.r.t. absolute compartment values), ‘least_squares_diff’ (least squares fit w.r.t. time-differences of compartment values)
  • alternative_guess (np.array, optional) – Alternative initial quess, different form the mean of the prior. Array in the same format as ‘flat_params’ in the result dictionary of a previous optimisation run.
  • use_mode_as_guess (bool, optional) – Initialise optimisation with mode instead of mean of the prior. Makes a difference for lognormal distributions.
  • tmp_file (optional, string) – If specified, name of a file to store the temporary best estimate of the global optimiser (as backup or for inspection) as numpy array file
  • load_backup_file (optional, string) – If specified, name of a file to restore the the state of the global optimiser
Returns:

output_dict – A dictionary containing the following keys for users:

x0: np.array

MAP estimates for the initial conditions

params_dict: dict

dictionary for MAP estimates for model parameters

control_params_dict: dict

dictionary for MAP estimates for control parameters (if requested)

-logp: float

Value of -logp at MAP.

Return type:

dict

Note

This function combines the functionality of latent_infer_parameters and latent_infer_control, which will be deprecated. To infer model parameters only, specify a fixed contactMatrix function. To infer control parameters only, specify a generator and do not specify priors for model parameters.

Examples

Here we list three examples, one for inferring all initial conditions along the fastest growing linear mode, one for inferring the initial conditions individually and a mixed one.

First, suppose we only observe Is out of (S, Ia, Is) and we wish to infer all compartmental values of S and Ia independently. For two age groups with population [2500, 7500],

>>> init_priors = {
        'independent':{
            'fltr': [True, True, True, True, False, False],
            'mean': [2400, 7400, 50, 50],
            'std': [200, 200, 200, 200],
            'bounds': [[2000, 2500], [7000, 7500], [0, 400], [0, 400]]
        }
    }

In the ‘fltr’ entry, we need a boolean array indicating which components of the full x0 = [S0[0], S0[1], Ia0[0], Ia0[1], Is0[0], Ia0[1]] array we are inferring. By setting fltr = [True, True, True, True, False, False], the inference algorithm will know that we are inferring all components of S0 and Ia0 but not Is0. Similar to inference for parameter values, we also assume a log-normal distribution for the priors for the initial conditions.

Next, if we are happy to assume that all our initial conditions lie along the fastest growing linear mode and we will only infer the coefficient of the mode, the init_priors dict would be,

>>> init_priors = {
        'lin_mode_coeff':{
            'fltr': [True, True, True, True, False, False],
            'mean': 100,
            'std': 100,
            'bounds': [1, 1000]
        }
    }

Note that the ‘fltr’ entry is still the same as before because we still only want to infer S and Ia, and the initial conditions for Is is fixed by the observation.

Finally, if we want to do a mixture of both (useful when some compartments have aligned with the fastest growing mode but others haven’t), we need to set the init_priors to be,

>>> init_priors = {
        'lin_mode_coeff': {
            'fltr': [True, True, False, False, False, False],
            'mean': 100,
            'std': 100,
            'bounds': [1, 1000]
        },
        'independent':{
            'fltr': [False, False, True, True, False, False],
            'mean': [50, 50],
            'std': [200, 200],
            'bounds': [0, 400], [0, 400]
        }
    }
latent_infer_control()

Compute the maximum a-posteriori (MAP) estimate of the change of control parameters for a SIR type model in lockdown with partially observed classes.

Parameters:latent_infer (see) –
Returns:output_dict – A dictionary containing the following keys for users:
map_params_dict: dict
dictionary for MAP estimates for control parameters
map_x0: np.array
MAP estimates for the initial conditions
-logp: float
Value of -logp at MAP.
Return type:dict

Note

This function just calls latent_infer (with the specified generator), will be deprecated.

latent_infer_mcmc()

Sample the posterior distribution of the initial conditions and all desired parameters, including control parameters, using ensemble MCMC. This requires the optional dependency emcee.

Parameters:
  • obs (2d numpy.array) – The observed trajectories with reduced number of variables (number of data points, (age groups * observed model classes))
  • fltr (2d numpy.array) – A matrix of shape (no. observed variables, no. total variables), such that obs_{ti} = fltr_{ij} * X_{tj}
  • Tf (float) – Total time of the trajectory
  • contactMatrix (callable) – A function that returns the contact matrix at time t (input).
  • param_priors (dict) – A dictionary for priors for the model parameters. See latent_infer_parameters for further explanations.
  • init_priors (dict) – A dictionary for priors for the initial conditions. See latent_infer_parameters for further explanations.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • verbose (bool, optional) – Set to True to see a progress bar for the sample generation. Default is False.
  • sampler (emcee.EnsembleSampler, optional) – Set to instance of the sampler (as returned by this function) to continue running the MCMC chains. Default is None (i.e. run a new chain).
  • nwalkers (int, optional) – The number of chains in the ensemble (should be at least 2*dim). Default is 2*dim.
  • walker_pos (np.array, optional) – The initial position of the walkers. If not specified, the function samples random positions from the prior.
  • nsamples (int, optional) – The number of samples per walker. Default is 1000.
  • nprocesses (int, optional) – The number of processes used to compute the likelihood for the walkers, needs pathos for values > 1. Default is the number of cpu cores if pathos is available, otherwise 1.
Returns:

sampler – This function returns the state of the sampler. To look at the chain of the internal flattened parameters, run sampler.get_chain(). Use this to judge whether the chain has sufficiently converged. Either rerun latent_infer_mcmc(…, sampler=sampler) to continue the chain or latent_infer_mcmc_process_result(…) to process the result.

Return type:

emcee.EnsembleSampler

Examples

For the structure of the model input parameters, in particular param_priors, init_priors, see the documentation of latent_infer. To start sampling the posterior, run for example

>>> sampler = estimator.latent_infer_mcmc(obs, fltr, Tf, param_priors, init_priors, contactMatrix=contactMatrix,
                                          verbose=True)

To judge the convergence of this chain, we can look at the trace plot of all the chains (for a moderate number of dimensions dim)

>>> fig, axes = plt.subplots(dim, sharex=True)
>>> samples = sampler.get_chain()
>>> for i in range(dim):
        ax = axes[i]
        ax.plot(samples[:, :, i], "k", alpha=0.3)
        ax.set_xlim(0, len(samples))
>>> axes[-1].set_xlabel("step number");

For more detailed convergence metrics, see the documentation of emcee. To continue running this chain, we can call this function again with the sampler as argument

>>> sampler = estimator.latent_infer_mcmc(obs, fltr, Tf, param_priors, init_priors, contactMatrix=contactMatrix,
                                          verbose=True, sampler=sampler)

This procudes 1000 additional samples in each chain. To process the results, call pyross.inference.latent_infer_mcmc_process_result.

latent_infer_mcmc_process_result()

Take the sampler generated by pyross.inference.latent_infer_mcmc and produce output dictionaries for further use in the pyross framework.

Parameters:
  • sampler (emcee.EnsembleSampler) – Output of mcmc_latent_inference.
  • obs (2d numpy.array) –
  • fltr (2d numpy.array) –
  • param_priors (dict) –
  • init_priors (dict) –
  • contactMatrix (callable, optional) –
  • generator (pyross.contactMatrix, optional) –
  • intervention_fun (callable, optional) –
  • flat (bool, optional) – This decides whether to return the samples as for each chain separately (False) or as as a combined list (True). Default is True.
  • discard (int, optional) – The number of initial samples to discard in each chain (to account for burn-in). Default is 0.
  • thin (int, optional) – Thin out the chain by taking only the n-tn element in each chain. Default is 1 (no thinning).
  • **catchall_kwargs (dict) – Catches further provided arguments and ignores them.
Returns:

output_samples – The processed posterior samples.

Return type:

list of dict (if flat=True), or list of list of dict (if flat=False)

latent_infer_nested_sampling()

Compute the log-evidence and weighted samples for the initial conditions and all desired parameters, including control parameters, for a SIR type model with partially observed classes. This function uses nested sampling as implemented in the dynesty Python package.

Parameters:
  • obs (2d numpy.array) – The observed trajectories with reduced number of variables (number of data points, (age groups * observed model classes))
  • fltr (2d numpy.array) – A matrix of shape (no. observed variables, no. total variables), such that obs_{ti} = fltr_{ij} * X_{tj}
  • Tf (float) – Total time of the trajectory
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input).
  • param_priors (dict) – A dictionary for priors for the model parameters. See latent_infer for further explanations.
  • init_priors (dict) – A dictionary for priors for the initial conditions. See latent_infer for further explanations.
  • contactMatrix – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to do inference in tangent space (might be less robust but a lot faster). Default is False.
  • verbose (bool, optional) – Set to True to see intermediate outputs from the nested sampling procedure.
  • nprocesses (int, optional) – The number of processes used for parallel evaluation of the likelihood.
  • queue_size (int, optional) – Size of internal queue of likelihood values, default is nprocesses if multiprocessing is used.
  • maxiter (int, optional) – The maximum number of iterations. Default is no limit.
  • maxcall (int, optional) – The maximum number of calls to the likelihood function. Default no limit.
  • dlogz (float, optional) – The iteration terminates if the estimated contribution of the remaining prior volume to the total evidence falls below this threshold. Default value is 1e-3 * (nlive - 1) + 0.01 if add_live==True, 0.01 otherwise.
  • n_effective (float, optional) – The iteration terminates if the number of effective posterior samples reaches this values. Default is no limit.
  • add_live (bool, optional) – Determines whether to add the remaining set of live points to the set of samples. Default is True.
  • sampler (dynesty.NestedSampler, optional) – Continue running an instance of a nested sampler until the termination criteria are met.
  • **dynesty_args

    Arguments passed through to the construction of the dynesty.NestedSampler constructor. Relevant entries are (this is not comprehensive, for details see the documentation of dynesty):

    nlive: int, optional
    The number of live points. Default is 500.
    bound: {‘none’, ‘single’, ‘multi’, ‘balls’, ‘cubes’}, optional
    Method used to approximately bound the prior using the current set of live points. Default is ‘multi’.
    sample: {‘auto’, ‘unif’, ‘rwalk’, ‘rstagger’, ‘slice’, ‘rslice’, ‘hslice’, callable}, optional
    Method used to sample uniformly within the likelihood constraint, conditioned on the provided bounds.
Returns:

sampler – The state of the sampler after termination of the nested sampling run.

Return type:

dynesty.NestedSampler

latent_infer_nested_sampling_process_result()

Take the sampler generated by pyross.inference.latent_infer_nested_sampling and produce output dictionaries for further use in the pyross framework. See there for additional description of parameters.

Parameters:
  • sampler (dynesty.NestedSampler) – Output of pyross.inference.latent_infer_nested_sampling.
  • obs (2d numpy.array) –
  • fltr (2d numpy.array) –
  • param_priors (dict) –
  • init_priors (dict) –
  • contactMatrix (callable, optional) –
  • generator (pyross.contactMatrix, optional) –
  • intervention_fun (callable, optional) –
  • **catchall_kwargs (dict) – Catches further provided arguments and ignores them.
Returns:

  • result (dynesty.Result) – The result of the nested sampling iteration. Relevant entries include:

    result.logz: list

    The progression of log-evidence estimates, use result.logz[-1] for the final estimate.

  • output_samples (list) – The processed weighted posterior samples.

latent_infer_parameters()

Compute the maximum a-posteriori (MAP) estimate of the parameters and the initial conditions of a SIR type model when the classes are only partially observed. Unobserved classes are treated as latent variables.

Parameters:latent_infer (see) –
Returns:output – Contains the following keys for users:
map_params_dict: dict
A dictionary for the MAP estimates for parameter values. The keys are the names of the parameters.
map_x0: np.array
The MAP estimate for the initial conditions.
-logp: float
The value of -logp at MAP.
Return type:dict

Note

This function just calls latent_infer (with fixed contactMatrix), will be deprecated.

latent_param_slice()

Samples the posterior and prior along a one-dimensional slice of the parameter space

Parameters:
  • obs (np.array) – The partially observed trajectory.
  • fltr (2d np.array) – The filter for the observation such that \(F_{ij} x_j (t) = obs_i(t)\)
  • Tf (float) – The total time of the trajectory.
  • infer_result (dict) – Dictionary returned by latent_infer
  • pos (np.array) – Position in parameter space around which the parameter slice is computed
  • direction (np.array) – Direction in parameter space in which the parameter slice is computed
  • scale (np.array) – Values by which the direction vector is scaled. Points evaluated are pos + scale * direction
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • inter_steps (int, optional) – Intermediate steps between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration.
  • nprocesses (int, optional) – The number of processes used to compute the likelihood for the walkers, needs pathos. Default is the number of cpu cores if pathos is available, otherwise 1.
Returns:

  • posterior (np.array) – posterior evaluated along the 1d slice
  • prior (np.array) – prior evaluated along the 1d slice

mcmc_inference()

Sample the posterior distribution of the epidemiological parameters using ensemble MCMC.

Note

This function has been replaced by pyross.inference.infer_mcmc and will be deleted in a future version of pyross. See there for a documentation of the function parameters.

mcmc_inference_process_result()

Take the sampler generated by mcmc_inference and produce output dictionaries for further use in the pyross framework.

Note

This function has been replaced by pyross.inference.infer_mcmc_process_result and will be deleted in a future version of pyross. See there for a documentation of the function parameters.

mcmc_latent_inference()

Sample the posterior distribution of the epidemiological parameters using ensemble MCMC.

Note

This function has been replaced by pyross.inference.latent_infer_mcmc and will be deleted in a future version of pyross. See there for a documentation of the function parameters.

mcmc_latent_inference_process_result()

Take the sampler generated by mcmc_latent_inference and produce output dictionaries for further use in the pyross framework.

Note

This function has been replaced by pyross.inference.latent_infer_mcmc_process_results and will be deleted in a future version of pyross. See there for a documentation of the function parameters.

minus_logp_red()

Computes -logp for a latent trajectory

Parameters:
  • parameters (dict) – A dictionary of parameter values, same as the ones required for initialisation.
  • x0 (numpy.array) – Initial conditions
  • obs (numpy.array) – The observed trajectory without the initial datapoint
  • fltr (boolean sequence or array) – True for observed and False for unobserved. e.g. if only Is is known for SIR with one age group, fltr = [False, False, True]
  • Tf (float) – The total time of the trajectory
  • contactMatrix (callable) – A function that returns the contact matrix at time t (input).
  • tangent (bool, optional) – Set to True to do inference in tangent space (might be less robust but a lot faster). Default is False.
Returns:

minus_logp – -log(p) for the observed trajectory with the given parameters and initial conditions

Return type:

float

nested_sampling_inference()

Run nested sampling for model parameters without latent variables.

Note

This function has been replaced by pyross.inference.infer_nested_sampling and will be deleted in a future version of pyross. See there for a documentation of the function parameters.

nested_sampling_inference_process_result()

Take the sampler generated by nested_sampling_inference and produce output dictionaries for further use in the pyross framework.

Note

This function has been replaced by pyross.inference.infer_nested_sampling_process_result and will be deleted in a future version of pyross. See there for a documentation of the function parameters.

nested_sampling_latent_inference()

Compute the log-evidence and weighted samples of the a-posteriori distribution of the parameters of a SIR type model with latent variables using nested sampling as implemented in the dynesty Python package.

Note

This function has been replaced by pyross.inference.latent_infer_nested_sampling and will be deleted in a future version of pyross. See there for a documentation of the function parameters.

nested_sampling_latent_inference_process_result()

Take the sampler generated by nested_sampling_latent_inference and produce output dictionaries for further use in the pyross framework.

Note

This function has been replaced by pyross.inference.latent_infer_nested_sampling_process_result and will be deleted in a future version of pyross. See there for a documentation of the function parameters.

obtain_minus_log_p()

Computes -logp of a full trajectory :param parameters: A dictionary for the model parameters. :type parameters: dict :param x: The full trajectory. :type x: np.array :param Tf: The time duration of the trajectory. :type Tf: float :param contactMatrix: A function that takes time (t) as an argument and returns the contactMatrix :type contactMatrix: callable :param tangent: Set to True to use tangent space inference. :type tangent: bool, optional

Returns:minus_logp – Value of -logp
Return type:float
robustness()

Robustness analysis in a two-dimensional slice of the parameter space, revealing neutral spaces as in https://doi.org/10.1073/pnas.1015814108.

Parameters:
  • FIM (2d numpy.array) – Fisher Information matrix of a stochastic model
  • FIM_det (2d numpy.array) – Fisher information matrix of the corresponding deterministic model
  • infer_result (dict) – Dictionary returned by latent_infer
  • param_pos_1 (int) – Position of ‘parameter 1’ in map_dict[‘flat_map’] for x-axis
  • param_pos_2 (int) – Position of ‘parameter 2’ in map_dict[‘flat_map’] for y-axis
  • range_1 (float) – Symmetric interval around parameter 1 for which robustness will be analysed. Absolute interval: ‘parameter 1’ +/- range_1
  • range_2 (float) – Symmetric interval around parameter 2 for which robustness will be analysed. Absolute interval: ‘parameter 2’ +/- range_2
  • resolution_1 (int) – Resolution of the meshgrid in x direction.
  • resolution_2 (int) – Resolution of the meshgrid in y direction. Default is resolution_2=resolution_1.
Returns:

  • ff (2d numpy.array) – shape=resolution_1 x resolution_2, meshgrid for x-axis
  • ss (2d numpy.array) – shape=resolution_1 x resolution_2, meshgrid for y-axis
  • Z_sto (2d numpy.array) – shape=resolution_1 x resolution_2, expected quadratic coefficient in the Taylor expansion of the likelihood of the stochastic model
  • Z_det (2d numpy.array) – shape=resolution_1 x resolution_2, expected quadratic coefficient in the Taylor expansion of the likelihood of the deterministic model

Examples

>>> from matplotlib import pyplot as plt
>>> from matplotlib import cm
>>>
>>> # positions 0 and 1 of infer_result['flat_params'] correspond to a scale parameter for alpha, and beta, respectively.
>>> ff, ss, Z_sto, Z_det = estimator.robustness(FIM, FIM_det, map_dict, 0, 1, 0.5, 0.01, 20)
>>> cmap = plt.cm.PuBu_r
>>> levels=11
>>> colors='black'
>>>
>>> c = plt.contourf(ff, ss, Z_sto, cmap=cmap, levels=levels) # heat map for the stochastic coefficient
>>> plt.contour(ff, ss, Z_sto, colors='black', levels=levels, linewidths=0.25)
>>> plt.contour(ff, ss, Z_det, colors=colors, levels=levels) # contour plot for the deterministic model
>>> plt.plot(infer_result['flat_params'][0], infer_result['flat_params'][1], 'o',
            color="#A60628", markersize=6) # the MAP estimate
>>> plt.colorbar(c)
>>> plt.xlabel(r'$lpha$ scale', fontsize=20, labelpad=10)
>>> plt.ylabel(r'$eta$', fontsize=20, labelpad=10)
>>> plt.show()
sample_gaussian()

Sample N samples of the parameters from the Gaussian centered at the MAP estimate with specified covariance cov.

Parameters:
  • N (int) – The number of samples.
  • map_estimate (dict) – The MAP estimate, e.g. as computed by inference.infer_parameters.
  • cov (np.array) – The covariance matrix of the flat parameters.
  • x (np.array) – The full trajectory.
  • Tf (float) – The total time of the trajectory.
  • contactMatrix (callable) – A function that returns the contact matrix at time t (input).
  • prior_dict (dict) – A dictionary containing priors. See examples.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
Returns:

samples – N samples of the Gaussian distribution.

Return type:

list of dict

sample_gaussian_latent()

Sample N samples of the parameters from the Gaussian centered at the MAP estimate with specified covariance cov.

Parameters:
  • N (int) – The number of samples.
  • obs (np.array) – The partially observed trajectory.
  • fltr (2d np.array) – The filter for the observation such that \(F_{ij} x_j (t) = obs_i(t)\)
  • Tf (float) – The total time of the trajectory.
  • infer_result (dict) – Dictionary returned by latent_infer
  • invcov (np.array) – The inverse covariance matrix of the flat parameters.
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • allow_negative (bool, optional) – Allow negative values of the sample parameters. If False, samples with negative paramters values are discarded and additional samples are drawn until the specified number N of samples is reached. Default is False.
  • inter_steps (int, optional) – Intermediate steps between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration.
  • nprocesses (int, optional) – The number of processes used to compute the likelihood for the walkers, needs pathos. Default is the number of cpu cores if pathos is available, otherwise 1.
Returns:

  • samples (list of np.array’s) – N samples of the Gaussian distribution (flat parameters).
  • posterior (np.array) – posterior evaluated along the 1d slice
  • prior (np.array) – prior evaluated along the 1d slice

sample_latent()

Samples the posterior and prior

Parameters:
  • obs (np.array) – The partially observed trajectory.
  • fltr (2d np.array) – The filter for the observation such that \(F_{ij} x_j (t) = obs_i(t)\)
  • Tf (float) – The total time of the trajectory.
  • infer_result (dict) – Dictionary returned by latent_infer
  • flat_params_list (list of np.array's) – Parameters for which the prior and posterior are sampled
  • contactMatrix (callable, optional) – A function that returns the contact matrix at time t (input). If specified, control parameters are not inferred. Either a contactMatrix or a generator must be specified.
  • generator (pyross.contactMatrix, optional) – A pyross.contactMatrix object that generates a contact matrix function with specified lockdown parameters. Either a contactMatrix or a generator must be specified.
  • intervention_fun (callable, optional) – The calling signature is intervention_func(t, **kwargs), where t is time and kwargs are other keyword arguments for the function. The function must return (aW, aS, aO), where aW, aS and aO are (2, M) arrays. The contact matrices are then rescaled as \(aW[0]_i CW_{ij} aW[1]_j\) etc. If not set, assume intervention that’s constant in time. See contactMatrix.constant_contactMatrix for details on the keyword parameters.
  • tangent (bool, optional) – Set to True to use tangent space inference. Default is False.
  • inter_steps (int, optional) – Intermediate steps between observations for the deterministic forward Euler integration. A higher number of intermediate steps will improve the accuracy of the result, but will make computations slower. Setting inter_steps=0 will fall back to the method accessible via det_method for the deterministic integration.
  • nprocesses (int, optional) – The number of processes used to compute the likelihood for the walkers, needs pathos. Default is the number of cpu cores if pathos is available, otherwise 1.
Returns:

  • posterior (np.array) – posterior evaluated along the 1d slice
  • prior (np.array) – prior evaluated along the 1d slice

sensitivity()

Computes sensitivity measures (not normalised) for 1) each individual parameter: from the diagonal elements of the FIM 2) incorporating parametric interactions: from the standard deviations derived from the FIM More on these interpretations can be found here: https://doi.org/10.1529/biophysj.104.053405 A larger entry translates into greater anticipated model sensitivity to changes in the parameter of interest.

Parameters:FIM (2d numpy.array) – The Fisher Information Matrix
Returns:
  • sensitivity_individual (numpy.array) – Sensitivity measure for individual parameters.
  • sensitivity_correlated (numpy.array) – Sensitivity measure incorporating parametric interactions.
set_contact_matrix()

Sets the internal contact matrix

Parameters:contactMatrix (callable) – A function that returns the contact matrix given time, with call signature contactMatrix(t).
set_det_method()

Sets the method used for deterministic integration for the SIR_type model

Parameters:
  • det_method (str) – The name of the integration method. Choose between ‘LSODA’ and ‘RK45’.
  • rtol (double, optional) – relative tolerance of the integrator (default 1e-3)
  • max_steps (int, optional) – Maximum number of integration steps (total) for the integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps are disregarded by the optimiser.
set_det_model()

Sets the internal deterministic model with given epidemiological parameters

Parameters:parameters (dict) – A dictionary of parameter values, same as the ones required for initialisation.
set_lyapunov_method()

Sets the method used for deterministic integration for the SIR_type model

Parameters:
  • lyapunov_method (str) – The name of the integration method. Choose between ‘LSODA’, ‘RK45’, ‘RK2’, ‘RK4’, and ‘euler’.
  • rtol (double, optional) – relative tolerance of the integrator (default 1e-3)
  • max_steps (int) – Maximum number of integration steps (total) for the integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps are disregarded by the optimiser.
set_params()

Sets epidemiological parameters used for evaluating -log(p)

Parameters:parameters (dict) – A dictionary containing all epidemiological parameters. Same keys as the one used to initialise the class.

Notes

Can use fill_params_dict to generate the full dictionary if only a few parameters are changed

Model

class pyross.inference.Model

Generic user-defined epidemic model.

To initialise the Model,

Parameters:
  • model_spec (dict) – A dictionary specifying the model. See Examples.
  • parameters (dict) – A dictionary containing the model parameters. All parameters can be float if not age-dependent, and np.array(M,) if age-dependent
  • M (int) – Number of age groups.
  • fi (np.array(M) or list) – Fraction of each age group.
  • Omega (int) – Total population.
  • steps (int, optional) – The number of internal integration steps performed between the observed points (not used in tangent space inference). For robustness, set steps to be large, lyapunov_method=’LSODA’. For speed, set steps to be small (~4), lyapunov_method=’euler’. For a combination of the two, choose something in between.
  • det_method (str, optional) – The integration method used for deterministic integration. Choose one of ‘LSODA’ and ‘RK45’. Default is ‘LSODA’.
  • lyapunov_method (str, optional) – The integration method used for the integration of the Lyapunov equation for the covariance. Choose one of ‘LSODA’, ‘RK45’, ‘RK2’, ‘RK4’ and ‘euler’. Default is ‘LSODA’.
  • rtol_det (float, optional) – relative tolerance for the deterministic integrator (default 1e-3)
  • rtol_lyapunov (float, optional) – relative tolerance for the Lyapunov-type integrator (default 1e-3)
  • max_steps_det (int, optional) – Maximum number of integration steps (total) for the deterministic integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_det are disregarded by the optimiser.
  • max_steps_lyapunov (int, optional) – Maximum number of integration steps (total) for the Lyapunov-type integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_lyapunov are disregarded by the optimiser.
  • parameter_mapping (python function, optional) – A user-defined function that maps the dictionary the parameters used for inference to a dictionary of parameters used in model_spec. Default is an identical mapping.
  • time_dep_param_mapping (python function, optional) – As parameter_mapping, but time-dependent. The user-defined function takes time as a second argument.
  • SIR_type for a table of all the methods (See) –

Examples

An example of model_spec and parameters for SIR class with a constant influx

>>> model_spec = {
        "classes" : ["S", "I"],
        "S" : {
            "constant"  : [ ["k"] ],
            "infection" : [ ["I", "S", "-beta"] ]
        },
        "I" : {
            "linear"    : [ ["I", "-gamma"] ],
            "infection" : [ ["I", "S", "beta"] ]
        }
    }
>>> parameters = {
        'beta': 0.1,
        'gamma': 0.1,
        'k': 1,
    }

Spp

class pyross.inference.Spp

This is a slightly more specific version of the class Model.

Spp is still supported for backward compatibility.

Model class is recommended over Spp for new users.

The Spp class works like Model but infection terms use a single class S

Parameters:
  • model_spec (dict) – A dictionary specifying the model. See Examples.
  • parameters (dict) – A dictionary containing the model parameters. All parameters can be float if not age-dependent, and np.array(M,) if age-dependent
  • M (int) – Number of age groups.
  • fi (np.array(M) or list) – Fraction of each age group.
  • Omega (int) – Total population.
  • steps (int, optional) – The number of internal integration steps performed between the observed points (not used in tangent space inference). For robustness, set steps to be large, lyapunov_method=’LSODA’. For speed, set steps to be small (~4), lyapunov_method=’euler’. For a combination of the two, choose something in between.
  • det_method (str, optional) – The integration method used for deterministic integration. Choose one of ‘LSODA’ and ‘RK45’. Default is ‘LSODA’.
  • lyapunov_method (str, optional) – The integration method used for the integration of the Lyapunov equation for the covariance. Choose one of ‘LSODA’, ‘RK45’, ‘RK2’ and ‘euler’. Default is ‘LSODA’.
  • rtol_det (float, optional) – relative tolerance for the deterministic integrator (default 1e-3)
  • rtol_lyapunov (float, optional) – relative tolerance for the Lyapunov-type integrator (default 1e-3)
  • max_steps_det (int, optional) – Maximum number of integration steps (total) for the deterministic integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_det are disregarded by the optimiser.
  • max_steps_lyapunov (int, optional) – Maximum number of integration steps (total) for the Lyapunov-type integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_lyapunov are disregarded by the optimiser.
  • parameter_mapping (python function, optional) – A user-defined function that maps the dictionary the parameters used for inference to a dictionary of parameters used in model_spec. Default is an identical mapping.
  • time_dep_param_mapping (python function, optional) – As parameter_mapping, but time-dependent. The user-defined function takes time as a second argument.
  • SIR_type for a table of all the methods (See) –

Examples

An example of model_spec and parameters for SIR class with a constant influx

>>> model_spec = {
        "classes" : ["S", "I"],
        "S" : {
            "constant"  : [ ["k"] ],
            "infection" : [ ["I", "-beta"] ]
        },
        "I" : {
            "linear"    : [ ["I", "-gamma"] ],
            "infection" : [ ["I", "beta"] ]
        }
    }
>>> parameters = {
        'beta': 0.1,
        'gamma': 0.1,
        'k': 1,
    }

SIR

class pyross.inference.SIR

Susceptible, Infected, Removed (SIR)

  • Ia: asymptomatic
  • Is: symptomatic

To initialise the SIR class,

Parameters:
  • parameters (dict) –

    Contains the following keys:

    alpha: float
    Ratio of asymptomatic carriers
    beta: float
    Infection rate upon contact
    gIa: float
    Recovery rate for asymptomatic
    gIs: float
    Recovery rate for symptomatic
    fsa: float
    Fraction by which symptomatic individuals do not self-isolate.
  • M (int) – Number of age groups
  • fi (float numpy.array) – Number of people in each age group divided by Omega.
  • Omega (float, optional) – System size parameter, e.g. total population. Default to 1.
  • steps (int, optional) – The number of internal integration steps performed between the observed points (not used in tangent space inference). For robustness, set steps to be large, lyapunov_method=’LSODA’. For speed, set steps to be small (~4), lyapunov_method=’euler’. For a combination of the two, choose something in between.
  • det_method (str, optional) – The integration method used for deterministic integration. Choose one of ‘LSODA’ and ‘RK45’. Default is ‘LSODA’.
  • lyapunov_method (str, optional) – The integration method used for the integration of the Lyapunov equation for the covariance. Choose one of ‘LSODA’, ‘RK45’, ‘RK2’, ‘RK4’ and ‘euler’. Default is ‘LSODA’.
  • rtol_det (float, optional) – relative tolerance for the deterministic integrator (default 1e-4)
  • rtol_lyapunov (float, optional) – relative tolerance for the Lyapunov-type integrator (default 1e-3)
  • max_steps_det (int, optional) – Maximum number of integration steps (total) for the deterministic integrator. Default: unlimited (represented as 0). Parameters for which the integrator reaches max_steps_det are disregarded by the optimiser.
  • max_steps_lyapunov (int, optional) – Maximum number of integration steps (total) for the Lyapunov-type integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_lyapunov are disregarded by the optimiser.

SEIR

class pyross.inference.SEIR

Susceptible, Exposed, Infected, Removed (SEIR)

  • Ia: asymptomatic
  • Is: symptomatic

To initialise the SEIR class,

Parameters:
  • parameters (dict) –

    Contains the following keys:

    alpha: float or np.array(M)
    Fraction of infected who are asymptomatic.
    beta: float
    Rate of spread of infection.
    gIa: float
    Rate of removal from asymptomatic individuals.
    gIs: float
    Rate of removal from symptomatic individuals.
    fsa: float
    Fraction by which symptomatic individuals do not self-isolate.
    gE: float
    rate of removal from exposed individuals.
  • M (int) – Number of age groups
  • fi (float numpy.array) – Number of people in each compartment divided by Omega
  • Omega (float, optional) – System size, e.g. total population. Default is 1.
  • steps (int, optional) – The number of internal integration steps performed between the observed points (not used in tangent space inference). For robustness, set steps to be large, lyapunov_method=’LSODA’. For speed, set steps to be small (~4), lyapunov_method=’euler’. For a combination of the two, choose something in between.
  • det_method (str, optional) – The integration method used for deterministic integration. Choose one of ‘LSODA’ and ‘RK45’. Default is ‘LSODA’.
  • lyapunov_method (str, optional) – The integration method used for the integration of the Lyapunov equation for the covariance. Choose one of ‘LSODA’, ‘RK45’, ‘RK2’, ‘RK4’ and ‘euler’. Default is ‘LSODA’.
  • rtol_det (float, optional) – relative tolerance for the deterministic integrator (default 1e-3)
  • rtol_lyapunov (float, optional) – relative tolerance for the Lyapunov-type integrator (default 1e-3)
  • max_steps_det (int, optional) – Maximum number of integration steps (total) for the deterministic integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_det are disregarded by the optimiser.
  • max_steps_lyapunov (int, optional) – Maximum number of integration steps (total) for the Lyapunov-type integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_lyapunov are disregarded by the optimiser.

SEAIRQ

class pyross.inference.SEAIRQ

Susceptible, Exposed, Asymptomatic and infected, Infected, Removed, Quarantined (SEAIRQ)

  • Ia: asymptomatic
  • Is: symptomatic
  • E: exposed
  • A: asymptomatic and infectious
  • Q: quarantined

To initialise the SEAIRQ class,

Parameters:
  • parameters (dict) –

    Contains the following keys:

    alpha: float or np.array(M)
    Fraction of infected who are asymptomatic.
    beta: float
    Rate of spread of infection.
    gIa: float
    Rate of removal from asymptomatic individuals.
    gIs: float
    Rate of removal from symptomatic individuals.
    gE: float
    rate of removal from exposed individuals.
    gA: float
    rate of removal from activated individuals.
    fsa: float
    Fraction by which symptomatic individuals do not self-isolate.
    tE: float
    testing rate and contact tracing of exposeds
    tA: float
    testing rate and contact tracing of activateds
    tIa: float
    testing rate and contact tracing of asymptomatics
    tIs: float
    testing rate and contact tracing of symptomatics
  • M (int) – Number of compartments
  • fi (float numpy.array) – Number of people in each compartment divided by Omega.
  • Omega (float, optional) – System size, e.g. total population. Default is 1.
  • steps (int, optional) – The number of internal integration steps performed between the observed points (not used in tangent space inference). For robustness, set steps to be large, lyapunov_method=’LSODA’. For speed, set steps to be small (~4), lyapunov_method=’euler’. For a combination of the two, choose something in between.
  • det_method (str, optional) – The integration method used for deterministic integration. Choose one of ‘LSODA’ and ‘RK45’. Default is ‘LSODA’.
  • lyapunov_method (str, optional) – The integration method used for the integration of the Lyapunov equation for the covariance. Choose one of ‘LSODA’, ‘RK45’, ‘RK2’, ‘RK4’ and ‘euler’. Default is ‘LSODA’.
  • rtol_det (float, optional) – relative tolerance for the deterministic integrator (default 1e-3)
  • rtol_lyapunov (float, optional) – relative tolerance for the Lyapunov-type integrator (default 1e-3)
  • max_steps_det (int, optional) – Maximum number of integration steps (total) for the deterministic integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_det are disregarded by the optimiser.
  • max_steps_lyapunov (int, optional) – Maximum number of integration steps (total) for the Lyapunov-type integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_lyapunov are disregarded by the optimiser.

SppQ

class pyross.inference.SppQ

User-defined epidemic model with quarantine stage.

This is a slightly more specific version of the class Model.

SppQ is still supported for backward compatibility.

Model class is recommended over SppQ for new users.

To initialise the SppQ model, …

Parameters:
  • model_spec (dict) – A dictionary specifying the model. See Examples.
  • parameters (dict) – A dictionary containing the model parameters. All parameters can be float if not age-dependent, and np.array(M,) if age-dependent
  • testRate (python function) – number of tests per day and age group
  • M (int) – Number of age groups.
  • fi (np.array(M) or list) – Fraction of each age group.
  • Omega (int) – Total population.
  • steps (int, optional) – The number of internal integration steps performed between the observed points (not used in tangent space inference). For robustness, set steps to be large, lyapunov_method=’LSODA’. For speed, set steps to be small (~4), lyapunov_method=’euler’. For a combination of the two, choose something in between.
  • det_method (str, optional) – The integration method used for deterministic integration. Choose one of ‘LSODA’ and ‘RK45’. Default is ‘LSODA’.
  • lyapunov_method (str, optional) – The integration method used for the integration of the Lyapunov equation for the covariance. Choose one of ‘LSODA’, ‘RK45’, ‘RK2’, ‘RK4’ and ‘euler’. Default is ‘LSODA’.
  • rtol_det (float, optional) – relative tolerance for the deterministic integrator (default 1e-3)
  • rtol_lyapunov (float, optional) – relative tolerance for the Lyapunov-type integrator (default 1e-3)
  • max_steps_det (int, optional) – Maximum number of integration steps (total) for the deterministic integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_det are disregarded by the optimiser.
  • max_steps_lyapunov (int, optional) – Maximum number of integration steps (total) for the Lyapunov-type integrator. Default: unlimited (represented as 0) Parameters for which the integrator reaches max_steps_lyapunov are disregarded by the optimiser.
  • parameter_mapping (python function, optional) – A user-defined function that maps the dictionary the parameters used for inference to a dictionary of parameters used in model_spec. Default is an identical mapping.
  • time_dep_param_mapping (python function, optional) – As parameter_mapping, but time-dependent. The user-defined function takes time as a second argument.
  • SIR_type for a table of all the methods (See) –

Examples

An example of model_spec and parameters for SIR class with random testing (without false positives/negatives) and quarantine

>>> model_spec = {
        "classes" : ["S", "I"],
        "S" : {
            "infection" : [ ["I", "-beta"] ]
        },
        "I" : {
            "linear"    : [ ["I", "-gamma"] ],
            "infection" : [ ["I", "beta"] ]
        },
        "test_pos"  : [ "p_falsepos", "p_truepos", "p_falsepos"] ,
        "test_freq" : [ "tf", "tf", "tf"]
    }
>>> parameters = {
        'beta': 0.1,
        'gamma': 0.1,
        'p_falsepos': 0
        'p_truepos': 1
        'tf': 1
    }