Quickstart#
pyeee
: a Python library for parameter screening of computational
models using Morris’ method of Elementary Effects and its extension of
Efficient or Sequential Elementary Effects by Cuntz, Mai et al.
(Water Res Research, 2015).
About pyeee#
pyeee
is a Python library for performing parameter screening of
computational models. It uses Morris’ method of Elementary Effects and
its extension, the so-called Efficient or Sequential Elementary Effects
published by
Cuntz, Mai et al. (2015) Computationally inexpensive identification of noninformative model parameters by sequential screening, Water Resources Research 51, 6417-6441, doi: 10.1002/2015WR016907.
pyeee
can be used with Python functions but also with external
programs, using for example the library partialwrap
. Function
evaluation can be distributed with Python’s multiprocessing
module or via the Message Passing Interface (MPI) mpi4py
.
Quick usage guide#
Simple Python function#
Consider the Ishigami-Homma function: \(y = sin(x_0) + a * sin(x_1)^2 + b * x_2^4 * sin(x_0)\).
Taking \(a = b = 1\) gives:
import numpy as np
def ishigami1(x):
return np.sin(x[0]) + np.sin(x[1])**2 + x[2]**4 * np.sin(x[0])
The three paramters \(x_0\), \(x_1\), \(x_2\) follow uniform distributions between \(-\pi\) and \(+\pi\).
Morris’ Elementary Effects can then be as:
from pyeee import screening
npars = 3
# lower boundaries
lb = np.ones(npars) * (-np.pi)
# upper boundaries
ub = np.ones(npars) * np.pi
# Elementary Effects
np.random.seed(seed=1023) # for reproducibility of examples
out = screening(ishigami1, lb, ub, 10) # mu*, mu, sigma
print("{:.1f} {:.1f} {:.1f}".format(*out[:, 0]))
# gives: 173.1 0.6 61.7
which gives the Elementary Effects mu*
.
Sequential Elementary Effects distinguish between informative and uninformative parameters using several times Morris’ Elementary Effects, returning a logical ndarray with True for the informative parameters and False for the uninformative parameters:
from pyeee import eee
# screen
np.random.seed(seed=1023) # for reproducibility of examples
out = eee(ishigami1, lb, ub, ntfirst=10)
print(out)
[ True False True]
Python function with extra parameters#
The function for for the routines in pyeee
must be of the form
\(func(x)\). Use Python’s partial()
from the
functools
module to pass other function parameters. For
example pass the parameters \(a\) and \(b\) to the
Ishigami-Homma function.
import numpy as np
from pyeee import eee
from functools import partial
def ishigami(x, a, b):
return np.sin(x[0]) + a * np.sin(x[1])**2 + b * x[2]**4 * np.sin(x[0])
def call_ishigami(func, a, b, x):
return func(x, a, b)
# Partialise function with fixed parameters
a = 0.5
b = 2.0
func = partial(call_ishigami, ishigami, a, b)
npars = 3
# lower boundaries
lb = np.ones(npars) * (-np.pi)
# upper boundaries
ub = np.ones(npars) * np.pi
# Elementary Effects
np.random.seed(seed=1023) # for reproducibility of examples
out = eee(func, lb, ub, ntfirst=10)
Figuratively speaking, partial()
passes a and
b to the function call_ishigami already during definition so that
eee
can then simply call it as func(x), where x is passed to
call_ishigami then as well.
Function wrappers#
We recommend to use our package partialwrap
for external
executables, which allows easy use of external programs and their
parallel execution. See the User Guide for
details. A trivial example is the use of partialwrap
for the
above function wrapping:
from partialwrap import function_wrapper
args = [a, b]
kwargs = {}
func = partial(func_wrapper, ishigami, args, kwargs)
# screen
out = eee(func, lb, ub, ntfirst=10)
Installation#
The easiest way to install is via pip:
pip install pyeee
or via conda:
conda install -c conda-forge pyeee
Requirements#
License#
pyeee
is distributed under the MIT License. See the
LICENSE file for details.
Copyright (c) 2019-2024 Matthias Cuntz, Juliane Mai
The project structure is based on a template provided by Sebastian Müller.