screening

Contents

screening#

Function for Morris’ method of Elementary Effects.

This function was written by Matthias Cuntz while at Institut National de Recherche pour l’Agriculture, l’Alimentation et l’Environnement (INRAE), Nancy, France.

copyright:

Copyright 2017-2022 Matthias Cuntz, see AUTHORS.rst for details.

license:

MIT License, see LICENSE for details.

The following functions are provided

ee(*args, **kwargs)

Wrapper function for screening().

screening(func, lb, ub, nt[, x0, mask, ...])

Parameter screening using Morris' method of Elementary Effects.

History
  • Written Dec 2017 by Matthias Cuntz (mc (at) macu (dot) de)

  • Output also (npara,3)-array for nt=1, Dec 2017, Matthias Cuntz

  • Removed call to external programs: use exe_wrappers, Jan 2018, Matthias Cuntz

  • Function can return multiple outputs, e.g. time series, Jan 2018, Matthias Cuntz

  • Python 3: map returns iterator, so list(map), Jul 2018, Fabio Gennaretti

  • Bug: default ntotal was not set if ntotal<0 (but nt instead), Dec 2019, Matthias Cuntz

  • Make numpy docstring format, Dec 2019, Matthias Cuntz

  • x0 optional, Jan 2020, Matthias Cuntz

  • Distinguish iterable and array_like parameter types; added seed keyword to screening/ee, Jan 2020, Matthias Cuntz

  • InputError does not exist, use TypeError, Feb 2020, Matthias Cuntz

  • Use new names of kwargs of morris_sampling and elementary_effects, Feb 2020, Matthias Cuntz

  • Make number of final trajectories an argument instead of a keyword argument, Feb 2020, Matthias Cuntz

  • Sample not only from uniform distribution but allow all distributions of scipy.stats, Mar 2020, Matthias Cuntz

  • Code refactoring, Sep 2021, Matthias Cuntz

  • More consistent docstrings, Jan 2022, Matthias Cuntz

  • Add functions with multiple outputs such as time series in Returns section, Apr 2024, Matthias Cuntz

ee(*args, **kwargs)[source]#

Wrapper function for screening().

screening(func, lb, ub, nt, x0=None, mask=None, nsteps=6, ntotal=None, dist=None, distparam=None, seed=None, processes=1, pool=None, verbose=0)[source]#

Parameter screening using Morris’ method of Elementary Effects.

Note, the input function must be callable as func(x).

Parameters:
  • func (callable) – Python function callable as func(x) with the function parameters x.

  • lb (array_like) – Lower bounds of parameters or lower fraction of percent point function ppf if distribution given. Be aware that the percent point function ppf of most continuous distributions is infinite at 0.

  • ub (array_like) – Upper bounds of parameters or upper fraction of percent point function ppf if distribution given. Be aware that the percent point function ppf of most continuous distributions is infinite at 1.

  • nt (int) – Number of trajectories used for screening.

  • x0 (array_like, optional) – Parameter values used with mask==0.

  • mask (array_like, optional) – Include (1,True) or exclude (0,False) parameters in screening (default: include all parameters).

  • nsteps (int, optional) – Number of steps along one trajectory (default: 6)

  • ntotal (int, optional) – Total number of sampled trajectories to select the nt most different trajectories. If None: max(nt**2,10*nt) (default: None)

  • dist (list, optional) – List of None or scipy.stats distribution objects for each factor having the method ppf, Percent Point Function (Inverse of CDF) (default: None). If None, the uniform distribution will be sampled from lower bound lb to upper bound ub. If dist is scipy.stats.uniform, the ppf will be sampled from the lower fraction given in lb and the upper fraction in ub. The sampling interval is then given by the parameters loc=lower and scale=interval=upper-lower in distparam. This means dist=None, lb=a, ub=b corresponds to lb=0, ub=1, dist=scipy.stats.uniform, distparam=[a,b-a]

  • distparam (list, optional) – List with tuples with parameters as required for dist (default: (0, 1)). All distributions of scipy.stats have location and scale parameters, at least. loc and scale are implemented as keyword arguments in scipy.stats. Other parameters such as the shape parameter of the gamma distribution must hence be given first, e.g. (shape,loc,scale) for the gamma distribution. distparam is ignored if dist is None. The percent point function ppf is called like this: dist(*distparam).ppf(x)

  • seed (int or array_like) – Seed for numpy’s random number generator (default: None).

  • processes (int, optional) – The number of processes to use to evaluate objective function and constraints (default: 1).

  • pool (schwimmbad pool object, optional) –

    Generic map function used from module schwimmbad, which provides, serial, multiprocessor, and MPI mapping functions (default: None). The pool will be chosen automatically if pool is None. It is chosen with:

    schwimmbad.choose_pool(mpi=True/False, processes=processes).
    

    MPI pools can only be opened and closed once. If you want to use screening several times in one program, then you have to choose the pool, pass it to screening, and later close the pool in the calling program.

  • verbose (int, optional) – Print progress report during execution if verbose>0 (default: 0).

Returns:

2D/3D-array with 3 entries per parameter and per function output, such as per time step in a time series:

  1. mean of absolute elementary effects over all nt trajectories (mu*)

  2. mean of elementary effects over all nt trajectories (mu)

  3. standard deviation of elementary effects over all nt trajectories (sigma) or zero if nt==1

Return type:

(nparameter, 3) ndarray or (nfuncout, nparameter, 3)

Examples

>>> from functools import partial
>>> import numpy as np
>>> from partialwrap import function_wrapper
>>> from pyjams.functions import fmorris
>>> seed = 1023
>>> np.random.seed(seed=seed)
>>> npars = 20
>>> beta0              = 0.
>>> beta1              = np.random.standard_normal(npars)
>>> beta1[:10]         = 20.
>>> beta2              = np.random.standard_normal((npars,npars))
>>> beta2[:6,:6]       = -15.
>>> beta3              = np.zeros((npars,npars,npars))
>>> beta3[:5,:5,:5]    = -10.
>>> beta4              = np.zeros((npars,npars,npars,npars))
>>> beta4[:4,:4,:4,:4] = 5.
>>> arg   = [beta0, beta1, beta2, beta3, beta4]
>>> kwarg = {}
>>> func  = partial(function_wrapper, fmorris, arg, kwarg)
>>> lb    = np.zeros(npars)
>>> ub    = np.ones(npars)
>>> nt      = 20
>>> ntotal  = 10*nt
>>> nsteps  = 6
>>> verbose = 0
>>> out = screening(func, lb, ub, nt, x0=None, mask=None,
...                 nsteps=nsteps, ntotal=ntotal, processes=4,
...                 verbose=verbose)
>>> print(out[0:3,0])
[60.7012889  67.33372626 48.46673528]