_sensitivity_analysis.py 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716
  1. import inspect
  2. from dataclasses import dataclass
  3. from typing import TYPE_CHECKING, Any
  4. from collections.abc import Callable
  5. import numpy as np
  6. from scipy.stats._common import ConfidenceInterval
  7. from scipy.stats._qmc import check_random_state
  8. from scipy.stats._resampling import BootstrapResult
  9. from scipy.stats import qmc, bootstrap
  10. from scipy._lib._array_api import xp_capabilities
  11. from scipy._lib._util import _transition_to_rng
  12. if TYPE_CHECKING:
  13. import numpy.typing as npt
  14. from scipy._lib._util import DecimalNumber, IntNumber
  15. __all__ = [
  16. 'sobol_indices'
  17. ]
  18. def f_ishigami(x: "npt.ArrayLike") -> "npt.NDArray[np.inexact[Any]]":
  19. r"""Ishigami function.
  20. .. math::
  21. Y(\mathbf{x}) = \sin x_1 + 7 \sin^2 x_2 + 0.1 x_3^4 \sin x_1
  22. with :math:`\mathbf{x} \in [-\pi, \pi]^3`.
  23. Parameters
  24. ----------
  25. x : array_like ([x1, x2, x3], n)
  26. Returns
  27. -------
  28. f : array_like (n,)
  29. Function evaluation.
  30. References
  31. ----------
  32. .. [1] Ishigami, T. and T. Homma. "An importance quantification technique
  33. in uncertainty analysis for computer models." IEEE,
  34. :doi:`10.1109/ISUMA.1990.151285`, 1990.
  35. """
  36. x = np.atleast_2d(x)
  37. f_eval = (
  38. np.sin(x[0])
  39. + 7 * np.sin(x[1])**2
  40. + 0.1 * (x[2]**4) * np.sin(x[0])
  41. )
  42. return f_eval
  43. def sample_A_B(
  44. n,
  45. dists,
  46. rng=None
  47. ):
  48. """Sample two matrices A and B.
  49. Uses a Sobol' sequence with 2`d` columns to have 2 uncorrelated matrices.
  50. This is more efficient than using 2 random draw of Sobol'.
  51. See sec. 5 from [1]_.
  52. Output shape is (d, n).
  53. References
  54. ----------
  55. .. [1] Saltelli, A., P. Annoni, I. Azzini, F. Campolongo, M. Ratto, and
  56. S. Tarantola. "Variance based sensitivity analysis of model
  57. output. Design and estimator for the total sensitivity index."
  58. Computer Physics Communications, 181(2):259-270,
  59. :doi:`10.1016/j.cpc.2009.09.018`, 2010.
  60. """
  61. d = len(dists)
  62. A_B = qmc.Sobol(d=2*d, seed=rng, bits=64).random(n).T
  63. A_B = A_B.reshape(2, d, -1)
  64. try:
  65. for d_, dist in enumerate(dists):
  66. A_B[:, d_] = dist.ppf(A_B[:, d_])
  67. except AttributeError as exc:
  68. message = "Each distribution in `dists` must have method `ppf`."
  69. raise ValueError(message) from exc
  70. return A_B
  71. def sample_AB(A: np.ndarray, B: np.ndarray) -> np.ndarray:
  72. """AB matrix.
  73. AB: rows of B into A. Shape (d, d, n).
  74. - Copy A into d "pages"
  75. - In the first page, replace 1st rows of A with 1st row of B.
  76. ...
  77. - In the dth page, replace dth row of A with dth row of B.
  78. - return the stack of pages
  79. """
  80. d, n = A.shape
  81. AB = np.tile(A, (d, 1, 1))
  82. i = np.arange(d)
  83. AB[i, i] = B[i]
  84. return AB
  85. def saltelli_2010(
  86. f_A: np.ndarray, f_B: np.ndarray, f_AB: np.ndarray
  87. ) -> tuple[np.ndarray, np.ndarray]:
  88. r"""Saltelli2010 formulation.
  89. .. math::
  90. S_i = \frac{1}{N} \sum_{j=1}^N
  91. f(\mathbf{B})_j (f(\mathbf{AB}^{(i)})_j - f(\mathbf{A})_j)
  92. .. math::
  93. S_{T_i} = \frac{1}{N} \sum_{j=1}^N
  94. (f(\mathbf{A})_j - f(\mathbf{AB}^{(i)})_j)^2
  95. Parameters
  96. ----------
  97. f_A, f_B : array_like (s, n)
  98. Function values at A and B, respectively
  99. f_AB : array_like (d, s, n)
  100. Function values at each of the AB pages
  101. Returns
  102. -------
  103. s, st : array_like (s, d)
  104. First order and total order Sobol' indices.
  105. References
  106. ----------
  107. .. [1] Saltelli, A., P. Annoni, I. Azzini, F. Campolongo, M. Ratto, and
  108. S. Tarantola. "Variance based sensitivity analysis of model
  109. output. Design and estimator for the total sensitivity index."
  110. Computer Physics Communications, 181(2):259-270,
  111. :doi:`10.1016/j.cpc.2009.09.018`, 2010.
  112. """
  113. # Empirical variance calculated using output from A and B which are
  114. # independent. Output of AB is not independent and cannot be used
  115. var = np.var([f_A, f_B], axis=(0, -1))
  116. # We divide by the variance to have a ratio of variance
  117. # this leads to eq. 2
  118. s = np.mean(f_B * (f_AB - f_A), axis=-1) / var # Table 2 (b)
  119. st = 0.5 * np.mean((f_A - f_AB) ** 2, axis=-1) / var # Table 2 (f)
  120. return s.T, st.T
  121. @dataclass
  122. class BootstrapSobolResult:
  123. first_order: BootstrapResult
  124. total_order: BootstrapResult
  125. @dataclass
  126. class SobolResult:
  127. first_order: np.ndarray
  128. total_order: np.ndarray
  129. _indices_method: Callable
  130. _f_A: np.ndarray
  131. _f_B: np.ndarray
  132. _f_AB: np.ndarray
  133. _A: np.ndarray | None = None
  134. _B: np.ndarray | None = None
  135. _AB: np.ndarray | None = None
  136. _bootstrap_result: BootstrapResult | None = None
  137. def bootstrap(
  138. self,
  139. confidence_level: "DecimalNumber" = 0.95,
  140. n_resamples: "IntNumber" = 999
  141. ) -> BootstrapSobolResult:
  142. """Bootstrap Sobol' indices to provide confidence intervals.
  143. Parameters
  144. ----------
  145. confidence_level : float, default: ``0.95``
  146. The confidence level of the confidence intervals.
  147. n_resamples : int, default: ``999``
  148. The number of resamples performed to form the bootstrap
  149. distribution of the indices.
  150. Returns
  151. -------
  152. res : BootstrapSobolResult
  153. Bootstrap result containing the confidence intervals and the
  154. bootstrap distribution of the indices.
  155. An object with attributes:
  156. first_order : BootstrapResult
  157. Bootstrap result of the first order indices.
  158. total_order : BootstrapResult
  159. Bootstrap result of the total order indices.
  160. See `BootstrapResult` for more details.
  161. """
  162. def statistic(idx):
  163. f_A_ = self._f_A[:, idx]
  164. f_B_ = self._f_B[:, idx]
  165. f_AB_ = self._f_AB[..., idx]
  166. return self._indices_method(f_A_, f_B_, f_AB_)
  167. n = self._f_A.shape[1]
  168. res = bootstrap(
  169. [np.arange(n)], statistic=statistic, method="BCa",
  170. n_resamples=n_resamples,
  171. confidence_level=confidence_level,
  172. bootstrap_result=self._bootstrap_result
  173. )
  174. self._bootstrap_result = res
  175. first_order = BootstrapResult(
  176. confidence_interval=ConfidenceInterval(
  177. res.confidence_interval.low[0], res.confidence_interval.high[0]
  178. ),
  179. bootstrap_distribution=res.bootstrap_distribution[0],
  180. standard_error=res.standard_error[0],
  181. )
  182. total_order = BootstrapResult(
  183. confidence_interval=ConfidenceInterval(
  184. res.confidence_interval.low[1], res.confidence_interval.high[1]
  185. ),
  186. bootstrap_distribution=res.bootstrap_distribution[1],
  187. standard_error=res.standard_error[1],
  188. )
  189. return BootstrapSobolResult(
  190. first_order=first_order, total_order=total_order
  191. )
  192. @xp_capabilities(np_only=True)
  193. @_transition_to_rng('random_state', replace_doc=False)
  194. def sobol_indices(
  195. *,
  196. func,
  197. n,
  198. dists=None,
  199. method='saltelli_2010',
  200. rng=None
  201. ):
  202. r"""Global sensitivity indices of Sobol'.
  203. Parameters
  204. ----------
  205. func : callable or dict(str, array_like)
  206. If `func` is a callable, function to compute the Sobol' indices from.
  207. Its signature must be::
  208. func(x: ArrayLike) -> ArrayLike
  209. with ``x`` of shape ``(d, n)`` and output of shape ``(s, n)`` where:
  210. - ``d`` is the input dimensionality of `func`
  211. (number of input variables),
  212. - ``s`` is the output dimensionality of `func`
  213. (number of output variables), and
  214. - ``n`` is the number of samples (see `n` below).
  215. Function evaluation values must be finite.
  216. If `func` is a dictionary, contains the function evaluations from three
  217. different arrays. Keys must be: ``f_A``, ``f_B`` and ``f_AB``.
  218. ``f_A`` and ``f_B`` should have a shape ``(s, n)`` and ``f_AB``
  219. should have a shape ``(d, s, n)``.
  220. This is an advanced feature and misuse can lead to wrong analysis.
  221. n : int
  222. Number of samples used to generate the matrices ``A`` and ``B``.
  223. Must be a power of 2. The total number of points at which `func` is
  224. evaluated will be ``n*(d+2)``.
  225. dists : list(distributions), optional
  226. List of each parameter's distribution. The distribution of parameters
  227. depends on the application and should be carefully chosen.
  228. Parameters are assumed to be independently distributed, meaning there
  229. is no constraint nor relationship between their values.
  230. Distributions must be an instance of a class with a ``ppf``
  231. method.
  232. Must be specified if `func` is a callable, and ignored otherwise.
  233. method : Callable or str, default: 'saltelli_2010'
  234. Method used to compute the first and total Sobol' indices.
  235. If a callable, its signature must be::
  236. func(f_A: np.ndarray, f_B: np.ndarray, f_AB: np.ndarray)
  237. -> Tuple[np.ndarray, np.ndarray]
  238. with ``f_A, f_B`` of shape ``(s, n)`` and ``f_AB`` of shape
  239. ``(d, s, n)``.
  240. These arrays contain the function evaluations from three different sets
  241. of samples.
  242. The output is a tuple of the first and total indices with
  243. shape ``(s, d)``.
  244. This is an advanced feature and misuse can lead to wrong analysis.
  245. rng : `numpy.random.Generator`, optional
  246. Pseudorandom number generator state. When `rng` is None, a new
  247. `numpy.random.Generator` is created using entropy from the
  248. operating system. Types other than `numpy.random.Generator` are
  249. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  250. .. versionchanged:: 1.15.0
  251. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  252. transition from use of `numpy.random.RandomState` to
  253. `numpy.random.Generator`, this keyword was changed from `random_state` to
  254. `rng`. For an interim period, both keywords will continue to work, although
  255. only one may be specified at a time. After the interim period, function
  256. calls using the `random_state` keyword will emit warnings. Following a
  257. deprecation period, the `random_state` keyword will be removed.
  258. Returns
  259. -------
  260. res : SobolResult
  261. An object with attributes:
  262. first_order : ndarray of shape (s, d)
  263. First order Sobol' indices.
  264. total_order : ndarray of shape (s, d)
  265. Total order Sobol' indices.
  266. And method:
  267. bootstrap(confidence_level: float, n_resamples: int)
  268. -> BootstrapSobolResult
  269. A method providing confidence intervals on the indices.
  270. See `scipy.stats.bootstrap` for more details.
  271. The bootstrapping is done on both first and total order indices,
  272. and they are available in `BootstrapSobolResult` as attributes
  273. ``first_order`` and ``total_order``.
  274. Notes
  275. -----
  276. The Sobol' method [1]_, [2]_ is a variance-based Sensitivity Analysis which
  277. obtains the contribution of each parameter to the variance of the
  278. quantities of interest (QoIs; i.e., the outputs of `func`).
  279. Respective contributions can be used to rank the parameters and
  280. also gauge the complexity of the model by computing the
  281. model's effective (or mean) dimension.
  282. .. note::
  283. Parameters are assumed to be independently distributed. Each
  284. parameter can still follow any distribution. In fact, the distribution
  285. is very important and should match the real distribution of the
  286. parameters.
  287. It uses a functional decomposition of the variance of the function to
  288. explore
  289. .. math::
  290. \mathbb{V}(Y) = \sum_{i}^{d} \mathbb{V}_i (Y) + \sum_{i<j}^{d}
  291. \mathbb{V}_{ij}(Y) + ... + \mathbb{V}_{1,2,...,d}(Y),
  292. introducing conditional variances:
  293. .. math::
  294. \mathbb{V}_i(Y) = \mathbb{\mathbb{V}}[\mathbb{E}(Y|x_i)]
  295. \qquad
  296. \mathbb{V}_{ij}(Y) = \mathbb{\mathbb{V}}[\mathbb{E}(Y|x_i x_j)]
  297. - \mathbb{V}_i(Y) - \mathbb{V}_j(Y),
  298. Sobol' indices are expressed as
  299. .. math::
  300. S_i = \frac{\mathbb{V}_i(Y)}{\mathbb{V}[Y]}
  301. \qquad
  302. S_{ij} =\frac{\mathbb{V}_{ij}(Y)}{\mathbb{V}[Y]}.
  303. :math:`S_{i}` corresponds to the first-order term which apprises the
  304. contribution of the i-th parameter, while :math:`S_{ij}` corresponds to the
  305. second-order term which informs about the contribution of interactions
  306. between the i-th and the j-th parameters. These equations can be
  307. generalized to compute higher order terms; however, they are expensive to
  308. compute and their interpretation is complex.
  309. This is why only first order indices are provided.
  310. Total order indices represent the global contribution of the parameters
  311. to the variance of the QoI and are defined as:
  312. .. math::
  313. S_{T_i} = S_i + \sum_j S_{ij} + \sum_{j,k} S_{ijk} + ...
  314. = 1 - \frac{\mathbb{V}[\mathbb{E}(Y|x_{\sim i})]}{\mathbb{V}[Y]}.
  315. First order indices sum to at most 1, while total order indices sum to at
  316. least 1. If there are no interactions, then first and total order indices
  317. are equal, and both first and total order indices sum to 1.
  318. .. warning::
  319. Negative Sobol' values are due to numerical errors. Increasing the
  320. number of points `n` should help.
  321. The number of sample required to have a good analysis increases with
  322. the dimensionality of the problem. e.g. for a 3 dimension problem,
  323. consider at minima ``n >= 2**12``. The more complex the model is,
  324. the more samples will be needed.
  325. Even for a purely additive model, the indices may not sum to 1 due
  326. to numerical noise.
  327. References
  328. ----------
  329. .. [1] Sobol, I. M.. "Sensitivity analysis for nonlinear mathematical
  330. models." Mathematical Modeling and Computational Experiment, 1:407-414,
  331. 1993.
  332. .. [2] Sobol, I. M. (2001). "Global sensitivity indices for nonlinear
  333. mathematical models and their Monte Carlo estimates." Mathematics
  334. and Computers in Simulation, 55(1-3):271-280,
  335. :doi:`10.1016/S0378-4754(00)00270-6`, 2001.
  336. .. [3] Saltelli, A. "Making best use of model evaluations to
  337. compute sensitivity indices." Computer Physics Communications,
  338. 145(2):280-297, :doi:`10.1016/S0010-4655(02)00280-1`, 2002.
  339. .. [4] Saltelli, A., M. Ratto, T. Andres, F. Campolongo, J. Cariboni,
  340. D. Gatelli, M. Saisana, and S. Tarantola. "Global Sensitivity Analysis.
  341. The Primer." 2007.
  342. .. [5] Saltelli, A., P. Annoni, I. Azzini, F. Campolongo, M. Ratto, and
  343. S. Tarantola. "Variance based sensitivity analysis of model
  344. output. Design and estimator for the total sensitivity index."
  345. Computer Physics Communications, 181(2):259-270,
  346. :doi:`10.1016/j.cpc.2009.09.018`, 2010.
  347. .. [6] Ishigami, T. and T. Homma. "An importance quantification technique
  348. in uncertainty analysis for computer models." IEEE,
  349. :doi:`10.1109/ISUMA.1990.151285`, 1990.
  350. Examples
  351. --------
  352. The following is an example with the Ishigami function [6]_
  353. .. math::
  354. Y(\mathbf{x}) = \sin x_1 + 7 \sin^2 x_2 + 0.1 x_3^4 \sin x_1,
  355. with :math:`\mathbf{x} \in [-\pi, \pi]^3`. This function exhibits strong
  356. non-linearity and non-monotonicity.
  357. Remember, Sobol' indices assumes that samples are independently
  358. distributed. In this case we use a uniform distribution on each marginals.
  359. >>> import numpy as np
  360. >>> from scipy.stats import sobol_indices, uniform
  361. >>> rng = np.random.default_rng()
  362. >>> def f_ishigami(x):
  363. ... f_eval = (
  364. ... np.sin(x[0])
  365. ... + 7 * np.sin(x[1])**2
  366. ... + 0.1 * (x[2]**4) * np.sin(x[0])
  367. ... )
  368. ... return f_eval
  369. >>> indices = sobol_indices(
  370. ... func=f_ishigami, n=1024,
  371. ... dists=[
  372. ... uniform(loc=-np.pi, scale=2*np.pi),
  373. ... uniform(loc=-np.pi, scale=2*np.pi),
  374. ... uniform(loc=-np.pi, scale=2*np.pi)
  375. ... ],
  376. ... rng=rng
  377. ... )
  378. >>> indices.first_order
  379. array([0.31637954, 0.43781162, 0.00318825])
  380. >>> indices.total_order
  381. array([0.56122127, 0.44287857, 0.24229595])
  382. Confidence interval can be obtained using bootstrapping.
  383. >>> boot = indices.bootstrap()
  384. Then, this information can be easily visualized.
  385. >>> import matplotlib.pyplot as plt
  386. >>> fig, axs = plt.subplots(1, 2, figsize=(9, 4))
  387. >>> _ = axs[0].errorbar(
  388. ... [1, 2, 3], indices.first_order, fmt='o',
  389. ... yerr=[
  390. ... indices.first_order - boot.first_order.confidence_interval.low,
  391. ... boot.first_order.confidence_interval.high - indices.first_order
  392. ... ],
  393. ... )
  394. >>> axs[0].set_ylabel("First order Sobol' indices")
  395. >>> axs[0].set_xlabel('Input parameters')
  396. >>> axs[0].set_xticks([1, 2, 3])
  397. >>> _ = axs[1].errorbar(
  398. ... [1, 2, 3], indices.total_order, fmt='o',
  399. ... yerr=[
  400. ... indices.total_order - boot.total_order.confidence_interval.low,
  401. ... boot.total_order.confidence_interval.high - indices.total_order
  402. ... ],
  403. ... )
  404. >>> axs[1].set_ylabel("Total order Sobol' indices")
  405. >>> axs[1].set_xlabel('Input parameters')
  406. >>> axs[1].set_xticks([1, 2, 3])
  407. >>> plt.tight_layout()
  408. >>> plt.show()
  409. .. note::
  410. By default, `scipy.stats.uniform` has support ``[0, 1]``.
  411. Using the parameters ``loc`` and ``scale``, one obtains the uniform
  412. distribution on ``[loc, loc + scale]``.
  413. This result is particularly interesting because the first order index
  414. :math:`S_{x_3} = 0` whereas its total order is :math:`S_{T_{x_3}} = 0.244`.
  415. This means that higher order interactions with :math:`x_3` are responsible
  416. for the difference. Almost 25% of the observed variance
  417. on the QoI is due to the correlations between :math:`x_3` and :math:`x_1`,
  418. although :math:`x_3` by itself has no impact on the QoI.
  419. The following gives a visual explanation of Sobol' indices on this
  420. function. Let's generate 1024 samples in :math:`[-\pi, \pi]^3` and
  421. calculate the value of the output.
  422. >>> from scipy.stats import qmc
  423. >>> n_dim = 3
  424. >>> p_labels = ['$x_1$', '$x_2$', '$x_3$']
  425. >>> sample = qmc.Sobol(d=n_dim, seed=rng).random(1024)
  426. >>> sample = qmc.scale(
  427. ... sample=sample,
  428. ... l_bounds=[-np.pi, -np.pi, -np.pi],
  429. ... u_bounds=[np.pi, np.pi, np.pi]
  430. ... )
  431. >>> output = f_ishigami(sample.T)
  432. Now we can do scatter plots of the output with respect to each parameter.
  433. This gives a visual way to understand how each parameter impacts the
  434. output of the function.
  435. >>> fig, ax = plt.subplots(1, n_dim, figsize=(12, 4))
  436. >>> for i in range(n_dim):
  437. ... xi = sample[:, i]
  438. ... ax[i].scatter(xi, output, marker='+')
  439. ... ax[i].set_xlabel(p_labels[i])
  440. >>> ax[0].set_ylabel('Y')
  441. >>> plt.tight_layout()
  442. >>> plt.show()
  443. Now Sobol' goes a step further:
  444. by conditioning the output value by given values of the parameter
  445. (black lines), the conditional output mean is computed. It corresponds to
  446. the term :math:`\mathbb{E}(Y|x_i)`. Taking the variance of this term gives
  447. the numerator of the Sobol' indices.
  448. >>> mini = np.min(output)
  449. >>> maxi = np.max(output)
  450. >>> n_bins = 10
  451. >>> bins = np.linspace(-np.pi, np.pi, num=n_bins, endpoint=False)
  452. >>> dx = bins[1] - bins[0]
  453. >>> fig, ax = plt.subplots(1, n_dim, figsize=(12, 4))
  454. >>> for i in range(n_dim):
  455. ... xi = sample[:, i]
  456. ... ax[i].scatter(xi, output, marker='+')
  457. ... ax[i].set_xlabel(p_labels[i])
  458. ... for bin_ in bins:
  459. ... idx = np.where((bin_ <= xi) & (xi <= bin_ + dx))
  460. ... xi_ = xi[idx]
  461. ... y_ = output[idx]
  462. ... ave_y_ = np.mean(y_)
  463. ... ax[i].plot([bin_ + dx/2] * 2, [mini, maxi], c='k')
  464. ... ax[i].scatter(bin_ + dx/2, ave_y_, c='r')
  465. >>> ax[0].set_ylabel('Y')
  466. >>> plt.tight_layout()
  467. >>> plt.show()
  468. Looking at :math:`x_3`, the variance
  469. of the mean is zero leading to :math:`S_{x_3} = 0`. But we can further
  470. observe that the variance of the output is not constant along the parameter
  471. values of :math:`x_3`. This heteroscedasticity is explained by higher order
  472. interactions. Moreover, an heteroscedasticity is also noticeable on
  473. :math:`x_1` leading to an interaction between :math:`x_3` and :math:`x_1`.
  474. On :math:`x_2`, the variance seems to be constant and thus null interaction
  475. with this parameter can be supposed.
  476. This case is fairly simple to analyse visually---although it is only a
  477. qualitative analysis. Nevertheless, when the number of input parameters
  478. increases such analysis becomes unrealistic as it would be difficult to
  479. conclude on high-order terms. Hence the benefit of using Sobol' indices.
  480. """
  481. rng = check_random_state(rng)
  482. n_ = int(n)
  483. if not (n_ & (n_ - 1) == 0) or n != n_:
  484. raise ValueError(
  485. "The balance properties of Sobol' points require 'n' "
  486. "to be a power of 2."
  487. )
  488. n = n_
  489. if not callable(method):
  490. indices_methods = {
  491. "saltelli_2010": saltelli_2010,
  492. }
  493. try:
  494. method = method.lower() # type: ignore[assignment]
  495. indices_method_ = indices_methods[method]
  496. except KeyError as exc:
  497. message = (
  498. f"{method!r} is not a valid 'method'. It must be one of"
  499. f" {set(indices_methods)!r} or a callable."
  500. )
  501. raise ValueError(message) from exc
  502. else:
  503. indices_method_ = method
  504. sig = inspect.signature(indices_method_)
  505. if set(sig.parameters) != {'f_A', 'f_B', 'f_AB'}:
  506. message = (
  507. "If 'method' is a callable, it must have the following"
  508. f" signature: {inspect.signature(saltelli_2010)}"
  509. )
  510. raise ValueError(message)
  511. def indices_method(f_A, f_B, f_AB):
  512. """Wrap indices method to ensure proper output dimension.
  513. 1D when single output, 2D otherwise.
  514. """
  515. return np.squeeze(indices_method_(f_A=f_A, f_B=f_B, f_AB=f_AB))
  516. if callable(func):
  517. if dists is None:
  518. raise ValueError(
  519. "'dists' must be defined when 'func' is a callable."
  520. )
  521. def wrapped_func(x):
  522. return np.atleast_2d(func(x))
  523. A, B = sample_A_B(n=n, dists=dists, rng=rng)
  524. AB = sample_AB(A=A, B=B)
  525. f_A = wrapped_func(A)
  526. if f_A.shape[1] != n:
  527. raise ValueError(
  528. "'func' output should have a shape ``(s, -1)`` with ``s`` "
  529. "the number of output."
  530. )
  531. def funcAB(AB):
  532. d, d, n = AB.shape
  533. AB = np.moveaxis(AB, 0, -1).reshape(d, n*d)
  534. f_AB = wrapped_func(AB)
  535. return np.moveaxis(f_AB.reshape((-1, n, d)), -1, 0)
  536. f_B = wrapped_func(B)
  537. f_AB = funcAB(AB)
  538. else:
  539. message = (
  540. "When 'func' is a dictionary, it must contain the following "
  541. "keys: 'f_A', 'f_B' and 'f_AB'."
  542. "'f_A' and 'f_B' should have a shape ``(s, n)`` and 'f_AB' "
  543. "should have a shape ``(d, s, n)``."
  544. )
  545. try:
  546. f_A, f_B, f_AB = map(lambda arr: arr.copy(), np.atleast_2d(
  547. func['f_A'], func['f_B'], func['f_AB']
  548. ))
  549. except KeyError as exc:
  550. raise ValueError(message) from exc
  551. if f_A.shape[1] != n or f_A.shape != f_B.shape or \
  552. f_AB.shape == f_A.shape or f_AB.shape[-1] % n != 0:
  553. raise ValueError(message)
  554. # Normalization by mean
  555. # Sobol', I. and Levitan, Y. L. (1999). On the use of variance reducing
  556. # multipliers in monte carlo computations of a global sensitivity index.
  557. # Computer Physics Communications, 117(1) :52-61.
  558. mean = np.mean([f_A, f_B], axis=(0, -1)).reshape(-1, 1)
  559. f_A -= mean
  560. f_B -= mean
  561. f_AB -= mean
  562. # Compute indices
  563. # Filter warnings for constant output as var = 0
  564. with np.errstate(divide='ignore', invalid='ignore'):
  565. first_order, total_order = indices_method(f_A=f_A, f_B=f_B, f_AB=f_AB)
  566. # null variance means null indices
  567. first_order[~np.isfinite(first_order)] = 0
  568. total_order[~np.isfinite(total_order)] = 0
  569. res = dict(
  570. first_order=first_order,
  571. total_order=total_order,
  572. _indices_method=indices_method,
  573. _f_A=f_A,
  574. _f_B=f_B,
  575. _f_AB=f_AB
  576. )
  577. if callable(func):
  578. res.update(
  579. dict(
  580. _A=A,
  581. _B=B,
  582. _AB=AB,
  583. )
  584. )
  585. return SobolResult(**res)