metacast.sensitivity_analyses
- Creation:
Author: Martin Grunnill Date: 2024-02-09
Description:
Submodules
Package Contents
Functions
|
Generate a Latin Hypercube sample, run model with sample and calculate PRCC for sampled parameters. |
|
Calculate Partial Correlation Coefficient PCC (default Rank 'Spearman' PRCC). |
|
Assess if a Latin-Hypercube (LH) sample size is reasonable. |
- metacast.sensitivity_analyses.lhs_prcc(parameters_df, sample_size, model_run_method, client=None, lhs_obj=None, other_samples_to_repeat=None, **kwargs)
Generate a Latin Hypercube sample, run model with sample and calculate PRCC for sampled parameters.
Latin Hypercube Sampling with Partial Rank Correlation Coefficients (PRCCS) can be used to evaluate sensitivity of a model to a parameter[1]. Note currently only supports uniform distribution.
- Parameters:
parameters_df (pandas.DataFrame or dictionary) – DataFrame outlining the boundaries for each parameter. Must contain fields ‘Lower Bound’ and ‘Upper Bound’. The name of the parameters is assumed to be in the index of the DataFrame. Alternatively, may be a dictionary that maps parameter names to percentage-point (inverse-CDF) functions.
sample_size (int) – Sample size of Latin Hypercube.
model_run_method (function) – Method of running model’s simulations. Must accept parameters as a single dictionary. Must output dictionary of input parameters and model results.
client (dask.distributed.Client (default None)) –
- Dask client for running simulations in parallel. If not given simulations are run serially with a tqdm progress
bar.
lhs_obj (scipy.stats.qmc.LatinHypercube, optional) – Pre-initialised Latin Hypercube sample generator. If not provided one is generated by within this function.
other_samples_to_repeat (pandas.DataFrame) – Samples to resampled and merged with LH samples.
kwargs – Key word arguments to be passed to model_run_method.
- Returns:
results_df (pandas.DataFrame) – Results of simulations preceded by parameters used in simulations.
prccs (pandas.DataFrame) – PRCC summary of the model’s sensitivity to its parameters.
References
- [1] Marino, S., Hogue, I. B., Ray, C. J., & Kirschner, D. E. (2008). A methodology for performing global uncertainty
and sensitivity analysis in systems biology. In Journal of Theoretical Biology (Vol. 254, Issue 1, pp. 178–196). https://doi.org/10.1016/j.jtbi.2008.04.011
- metacast.sensitivity_analyses.calculate_prcc(results_and_sample_df, parameter, output, covariables, method='spearman')
Calculate Partial Correlation Coefficient PCC (default Rank ‘Spearman’ PRCC).
A wrapper for pingouin.partial_corr. Partial Rank Correlation Coefficients (PRCCS) can be used to evaluate sensitivity of a model to a parameter[1].
- Parameters:
results_and_sample_df (pandas.DataFrame) – DataFrame of results and parameters for calculating PCC.
parameter (string) – Parameter for which PCC will be calculated.
output (string) – Model output for which PCC will be calculated.
covariables (list of strings) – Parameters whose effects will be discounted.
method (string, default 'spearman' (Rank correlations)) – Form of PCC see documentation of pingouin.partial_corr.
- Returns:
Partial Corelation Coefficient of parameter and output.
- Return type:
pandas.DataFrame
References
- [1] Marino, S., Hogue, I. B., Ray, C. J., & Kirschner, D. E. (2008). A methodology for performing global uncertainty
and sensitivity analysis in systems biology. In Journal of Theoretical Biology (Vol. 254, Issue 1, pp. 178–196). https://doi.org/10.1016/j.jtbi.2008.04.011
- metacast.sensitivity_analyses.determine_lh_sample_size(parameters_df, model_run_method, start_n, repeats_per_n, std_aim, LHS_PRCC_method, save_dir, attempts_to_make=float('inf'), n_increase_addition=None, n_increase_multiple=None, y0=None, other_samples_to_repeat=None, max_workers=None)
Assess if a Latin-Hypercube (LH) sample size is reasonable.
- Method is as follows:
For a given LH sample size draw and simulate a model through repeats_per_n number of LHSs.
Determine standard deviation in PRCCs. If any PRCC is greater than std_aim increase LH sample size by either addittion (n_increase_addition) or multiplication (n_increase_multiple) and return to 1.
- Parameters:
parameters_df (pd.DataFrame) – DataFrame outlining the boundaries for each parameter. Must contain fields ‘Lower Bound’ and ‘Upper Bound’. Name of the parameters is assumed to be in the index.
model_run_method (function) – Method for simulating model.
start_n (int) – Suggested sample size to start from.
repeats_per_n (int) – Number of repeats per sample size.
std_aim (float) – Standard deviation goal.
LHS_PRCC_method (function) – Method for drawing LHs simulating models through LHs and determining PRCC.
save_dir (string) – Directory for saving results.
attempts_to_make (int) – Number of attempts to make.
n_increase_addition (int, mutually exclusive with n_increase_multiple) – Increase in sample size, by addition, to be made if standard deviation in any paramters PRCC is greater than std_aim.
n_increase_multiple (int, mutually exclusive with n_increase_addition) – Increase in sample size, by multiplication, to be made if standard deviation in any paramters PRCC is greater than std_aim.
y0 (numpy.array, optional) – Intial values of state variables.
other_samples_to_repeat (pandas.DataFrame) – Samples to resampled and merged with LH samples.
max_workers (int, optional) – Number workers if LHS_PRCC method is parallelized.
- Returns:
Sample size assessed as being reasonable.
- Return type:
int