soboljansen {sensitivity} | R Documentation |
Monte Carlo Estimation of Sobol' Indices (improved formulas of Jansen (1999) and Saltelli et al. (2010))
Description
soboljansen
implements the Monte Carlo estimation of
the Sobol' indices for both first-order and total indices at the same
time (alltogether 2p
indices), at a total cost of (p+2)
\times n
model evaluations. These are called the Jansen estimators.
Usage
soboljansen(model = NULL, X1, X2, nboot = 0, conf = 0.95, ...)
## S3 method for class 'soboljansen'
tell(x, y = NULL, return.var = NULL, ...)
## S3 method for class 'soboljansen'
print(x, ...)
## S3 method for class 'soboljansen'
plot(x, ylim = c(0, 1), y_col = NULL, y_dim3 = NULL, ...)
## S3 method for class 'soboljansen'
plotMultOut(x, ylim = c(0, 1), ...)
## S3 method for class 'soboljansen'
ggplot(data, mapping = aes(), ylim = c(0, 1), y_col = NULL,
y_dim3 = NULL, ..., environment = parent.frame())
Arguments
model |
a function, or a model with a |
X1 |
the first random sample. |
X2 |
the second random sample. |
nboot |
the number of bootstrap replicates. |
conf |
the confidence level for bootstrap confidence intervals. |
x |
a list of class |
data |
a list of class |
y |
a vector of model responses. |
return.var |
a vector of character strings giving further
internal variables names to store in the output object |
ylim |
y-coordinate plotting limits. |
y_col |
an integer defining the index of the column of |
y_dim3 |
an integer defining the index in the third dimension of
|
mapping |
Default list of aesthetic mappings to use for plot. If not specified, must be supplied in each layer added to the plot. |
environment |
[Deprecated] Used prior to tidy evaluation. |
... |
for |
Details
This estimator is good for large first-order indices, and (large and small) total indices.
This version of soboljansen
also supports matrices and three-dimensional
arrays as output of model
. If the model output is a matrix or an array,
V
, S
and T
are matrices or arrays as well (depending on the
type of y
and the value of nboot
).
The bootstrap outputs V.boot
, S.boot
and T.boot
can only be
returned if the model output is a vector (using argument return.var
). For
matrix or array output, these objects can't be returned.
Value
soboljansen
returns a list of class "soboljansen"
, containing all
the input arguments detailed before, plus the following components:
call |
the matched call. |
X |
a |
y |
either a vector, a matrix or a three-dimensional array of model
responses (depends on the output of |
V |
the estimations of Variances of the Conditional Expectations
(VCE) with respect to each factor and also with respect to the
complementary set of each factor ("all but |
S |
the estimations of the Sobol' first-order indices. |
T |
the estimations of the Sobol' total sensitivity indices. |
Users can ask more ouput variables with the argument
return.var
(for example, bootstrap outputs V.boot
,
S.boot
and T.boot
).
Author(s)
Bertrand Iooss, with contributions from Frank Weber (2016)
References
M.J.W. Jansen, 1999, Analysis of variance designs for model output, Computer Physics Communication, 117, 35–43.
A. Saltelli, P. Annoni, I. Azzini, F. Campolongo, M. Ratto and S. Tarantola, 2010, Variance based sensitivity analysis of model output. Design and estimator for the total sensitivity index, Computer Physics Communications 181, 259–270.
See Also
sobol, sobol2002, sobolSalt, sobol2007, sobolmartinez, sobolEff, sobolmara,sobolMultOut
Examples
# Test case : the non-monotonic Sobol g-function
# The method of sobol requires 2 samples
# There are 8 factors, all following the uniform distribution
# on [0,1]
library(boot)
n <- 1000
X1 <- data.frame(matrix(runif(8 * n), nrow = n))
X2 <- data.frame(matrix(runif(8 * n), nrow = n))
# sensitivity analysis
x <- soboljansen(model = sobol.fun, X1, X2, nboot = 100)
print(x)
plot(x)
library(ggplot2)
ggplot(x)
# Only for demonstration purposes: a model function returning a matrix
sobol.fun_matrix <- function(X){
res_vector <- sobol.fun(X)
cbind(res_vector, 2 * res_vector)
}
x_matrix <- soboljansen(model = sobol.fun_matrix, X1, X2)
plot(x_matrix, y_col = 2)
title(main = "y_col = 2")
# Also only for demonstration purposes: a model function returning a
# three-dimensional array
sobol.fun_array <- function(X){
res_vector <- sobol.fun(X)
res_matrix <- cbind(res_vector, 2 * res_vector)
array(data = c(res_matrix, 5 * res_matrix),
dim = c(length(res_vector), 2, 2))
}
x_array <- soboljansen(model = sobol.fun_array, X1, X2)
plot(x_array, y_col = 2, y_dim3 = 2)
title(main = "y_col = 2, y_dim3 = 2")