AdaptiveDirectionalSampling

class AdaptiveDirectionalSampling(*args)

Adaptative directional simulation.

Parameters
eventRandomVector

Event we are computing the probability of.

rootStrategyRootStrategy, optional

Strategy adopted to evaluate the intersections of each direction with the limit state function and take into account the contribution of the direction to the event probability. Set to SafeAndSlow by default.

samplingStrategySamplingStrategy, optional

Strategy adopted to sample directions. Set to RandomDirection by default.

See also

EventSimulation

Notes

Let \mathcal D_f denote the failure domain defined as \mathcal D_f = \{\ux \in \mathbb R^{n_X} | g(\ux) \leq 0\}, where \ux are realization of the random vector \uX and g is the limit-state function as defined elsewhere in the documentation.

The purpose of the ADS-2 algorithm and its variants is to estimate the following probability:

P_f = \int_{\mathcal D_f} f_{\uX}(\ux)\di{\ux} \\
    = \int_{\mathbb R^{n_X}} \mathbf{1}_{\{g(\ux) \:\leq 0\: \}}\,f_{\uX}(\ux)\di{\ux} \\
    = \Prob{\{g(\uX) \leq 0\}}.

Principles

The ADS-2 method [munoz2011] combines the stratified and directional sampling concepts. Stratified sampling consists in splitting the support of the random vector \ux into m mutually exclusive and collectively exhaustive subsets. Here, ADS-2 splits the standard space into m = 2^d quadrants, where d is the dimension of the random vector \uX. Stratified sampling is often run in two steps: (i) a learning step is used for polling the input space and detect the subsets that contribute most to the probability and (ii) an estimation step is used for estimating the probability by weighted sampling (some subsets are more sampled than the others). Directional sampling uses the spheric symmetry of the standard space for estimating the failure probability as the average of conditional probabilities calculated on directions drawn at random in the standard space.

The learning step uses an a priori number of random directions that is uniformly distributed over the quadrants, meaning the weights are as follows:

\omega^1_i = \frac{1}{m}, \quad i = 1, \ldots, m.

Directional sampling is used for estimating the failure probability in each quadrant:

\hat P_i^{DS} = \Prob{\{g(\uX) \leq 0\} \mid \uX \in \mathbb{Q}_i},\,i = 1, \ldots, m.

and the corresponding estimation variances are denoted as \sigma_i^{DS\,2}. These probabilities are estimated using the same number N^0_i of random directions per quadrant as told by the uniform weights distribution.

The probability of interest is then computed as a weighted average of the previously defined conditional probabilities:

\hat P_f = \sum\limits_{i=1}^m \omega_i \hat P_i^{DS}

where \hat P_i^{DS} is the conditional probability estimator in the i-th quadrant. The corresponding variance of the stratified estimator reads:

\sigma^2 = \frac{1}{N_l} \left( \sum\limits_{i=1}^m \omega_i \sigma_i^{DS} \right) ^2

where \sigma_i^{DS\,2} is the variance of the conditional probability estimator in the i-th quadrant.

At the end of the learning step, the weights \omega_i are updated so as to minimize the stratified estimator variance. Indeed, it can be shown that the updated weights:

\omega^2_i = \frac{\omega^1_i \sigma_i^{DS}}{\sum\limits_{j=1}^m \omega^1_j \sigma_j^{DS}}, i = 1, \ldots, m,

minimize the final estimation variance in eqref{eq:pf_est_sda2_var}. Note that some weights might be zero (due to a somewhat arbitrary rounding of the conditional probabilities’ estimation variance). The quadrants associated with a zero-weight will not be sampled in the estimation step.

Eventually, the estimation step proceeds in essentially the same way as the learning step with different weights for the quadrants though. eqref{eq:pf_est_sda2} and eqref{eq:pf_est_sda2_var} are used for evaluating the final probability probability estimate and its variance.

The computational budget per step is parametrized by a fraction \gamma_l, l = 1,\,2 of the total budget N, such that \gamma_1 + \gamma_2 = 1. The number of directions sampled in quadrant i at step l is then defined as follows:

N^l_i = N * \gamma_l * \omega_i.

The number of evaluation of the limit-state function g is of course greater than the total budget N since directional sampling is used.

Variants

The ADS-2+ variant performs a dimension reduction step after the learning step for reducing the number of stratified quadrants. The statistic \tilde T_k aggregates the sensitivity of expectation along dimension k. It is defined as follows:

\tilde T_k = \sum\limits_{i_l \in \lbrace -1,1 \rbrace,l \neq k} \lvert \tilde I_{(i_1,\dots,i_{k-1},-1,i_{k+1},\dots,i_p)} - \tilde I_{(i_1,\dots,i_{k-1},1,i_{k+1},\dots,i_p)} \rvert.

It is used for ranking the contributions of the quadrants. Then, only the d' < d most influential variables according to \tilde T_k are stratified, leaving the remaining variables simulated without stratification. The corresponding quadrants will not be sampled.

The DP-ADS-2 variant combines the ADS method with a rotation of the quadrants. The idea is to get a possible design point (available e.g. after a preliminary FORM analysis) on the bisector of one of the quadrants to make the stratification even more efficient and thus save some evaluations of the model.

This 2-step algorithm can be generalized to L > 2 steps by adding more than one learning step. For now, only ADS-2 is implemented.

Methods

drawProbabilityConvergence(*args)

Draw the probability convergence at a given level.

getBlockSize()

Accessor to the block size.

getClassName()

Accessor to the object’s name.

getConvergenceStrategy()

Accessor to the convergence strategy.

getEvent()

Accessor to the event.

getGamma()

Gamma accessor.

getId()

Accessor to the object’s id.

getMaximumCoefficientOfVariation()

Accessor to the maximum coefficient of variation.

getMaximumOuterSampling()

Accessor to the maximum sample size.

getMaximumStandardDeviation()

Accessor to the maximum standard deviation.

getMaximumStratificationDimension()

Maximum stratification dimension accessor.

getName()

Accessor to the object’s name.

getPartialStratification()

Partial stratification accessor.

getQuadrantOrientation()

Quadrant orientation accessor.

getResult()

Accessor to the results.

getRootStrategy()

Get the root strategy.

getSamplingStrategy()

Get the direction sampling strategy.

getShadowedId()

Accessor to the object’s shadowed id.

getTStatistic()

T statistic accessor.

getVerbose()

Accessor to verbosity.

getVisibility()

Accessor to the object’s visibility state.

hasName()

Test if the object is named.

hasVisibleName()

Test if the object has a distinguishable name.

run()

Launch simulation.

setBlockSize(blockSize)

Accessor to the block size.

setConvergenceStrategy(convergenceStrategy)

Accessor to the convergence strategy.

setGamma(gamma)

Gamma accessor.

setMaximumCoefficientOfVariation(…)

Accessor to the maximum coefficient of variation.

setMaximumOuterSampling(maximumOuterSampling)

Accessor to the maximum sample size.

setMaximumStandardDeviation(…)

Accessor to the maximum standard deviation.

setMaximumStratificationDimension(…)

Maximum stratification dimension accessor.

setName(name)

Accessor to the object’s name.

setPartialStratification(partialStratification)

Partial stratification accessor.

setProgressCallback(*args)

Set up a progress callback.

setQuadrantOrientation(quadrantOrientation)

Quadrant orientation accessor.

setRootStrategy(rootStrategy)

Set the root strategy.

setSamplingStrategy(samplingStrategy)

Set the direction sampling strategy.

setShadowedId(id)

Accessor to the object’s shadowed id.

setStopCallback(*args)

Set up a stop callback.

setVerbose(verbose)

Accessor to verbosity.

setVisibility(visible)

Accessor to the object’s visibility state.

__init__(*args)

Initialize self. See help(type(self)) for accurate signature.

drawProbabilityConvergence(*args)

Draw the probability convergence at a given level.

Parameters
levelfloat, optional

The probability convergence is drawn at this given confidence length level. By default level is 0.95.

Returns
grapha Graph

probability convergence graph

getBlockSize()

Accessor to the block size.

Returns
blockSizeint

Number of terms in the probability simulation estimator grouped together. It is set by default to 1.

getClassName()

Accessor to the object’s name.

Returns
class_namestr

The object class name (object.__class__.__name__).

getConvergenceStrategy()

Accessor to the convergence strategy.

Returns
storage_strategyHistoryStrategy

Storage strategy used to store the values of the probability estimator and its variance during the simulation algorithm.

getEvent()

Accessor to the event.

Returns
eventRandomVector

Event we want to evaluate the probability.

getGamma()

Gamma accessor.

The computational budget per step \gamma_l.

Returns
gammaPoint

Gamma value.

getId()

Accessor to the object’s id.

Returns
idint

Internal unique identifier.

getMaximumCoefficientOfVariation()

Accessor to the maximum coefficient of variation.

Returns
coefficientfloat

Maximum coefficient of variation of the simulated sample.

getMaximumOuterSampling()

Accessor to the maximum sample size.

Returns
outerSamplingint

Maximum number of groups of terms in the probability simulation estimator.

getMaximumStandardDeviation()

Accessor to the maximum standard deviation.

Returns
sigmafloat, \sigma > 0

Maximum standard deviation of the estimator.

getMaximumStratificationDimension()

Maximum stratification dimension accessor.

Returns
maxint

Maximum stratification dimension.

getName()

Accessor to the object’s name.

Returns
namestr

The name of the object.

getPartialStratification()

Partial stratification accessor.

Returns
partialStratificationbool

Partial stratification.

getQuadrantOrientation()

Quadrant orientation accessor.

Returns
orientationPoint

Quadrant orientation.

getResult()

Accessor to the results.

Returns
resultsSimulationResult

Structure containing all the results obtained after simulation and created by the method run().

getRootStrategy()

Get the root strategy.

Returns
strategyRootStrategy

Root strategy adopted.

getSamplingStrategy()

Get the direction sampling strategy.

Returns
strategySamplingStrategy

Direction sampling strategy adopted.

getShadowedId()

Accessor to the object’s shadowed id.

Returns
idint

Internal unique identifier.

getTStatistic()

T statistic accessor.

The statistic \tilde T_k aggregates the sensitivity of expectation.

Returns
gammaPoint

T statistic value.

getVerbose()

Accessor to verbosity.

Returns
verbosity_enabledbool

If True, the computation is verbose. By default it is verbose.

getVisibility()

Accessor to the object’s visibility state.

Returns
visiblebool

Visibility flag.

hasName()

Test if the object is named.

Returns
hasNamebool

True if the name is not empty.

hasVisibleName()

Test if the object has a distinguishable name.

Returns
hasVisibleNamebool

True if the name is not empty and not the default one.

run()

Launch simulation.

Notes

It launches the simulation and creates a SimulationResult, structure containing all the results obtained after simulation. It computes the probability of occurence of the given event by computing the empirical mean of a sample of size at most outerSampling * blockSize, this sample being built by blocks of size blockSize. It allows to use efficiently the distribution of the computation as well as it allows to deal with a sample size > 2^{32} by a combination of blockSize and outerSampling.

setBlockSize(blockSize)

Accessor to the block size.

Parameters
blockSizeint, blockSize \geq 1

Number of terms in the probability simulation estimator grouped together. It is set by default to 1.

Notes

For Monte Carlo, LHS and Importance Sampling methods, this allows to save space while allowing multithreading, when available we recommend to use the number of available CPUs; for the Directional Sampling, we recommend to set it to 1.

setConvergenceStrategy(convergenceStrategy)

Accessor to the convergence strategy.

Parameters
storage_strategyHistoryStrategy

Storage strategy used to store the values of the probability estimator and its variance during the simulation algorithm.

setGamma(gamma)

Gamma accessor.

The computational budget per step \gamma_l.

Parameters
gammasequence of float

Gamma value.

setMaximumCoefficientOfVariation(maximumCoefficientOfVariation)

Accessor to the maximum coefficient of variation.

Parameters
coefficientfloat

Maximum coefficient of variation of the simulated sample.

setMaximumOuterSampling(maximumOuterSampling)

Accessor to the maximum sample size.

Parameters
outerSamplingint

Maximum number of groups of terms in the probability simulation estimator.

setMaximumStandardDeviation(maximumStandardDeviation)

Accessor to the maximum standard deviation.

Parameters
sigmafloat, \sigma > 0

Maximum standard deviation of the estimator.

setMaximumStratificationDimension(maximumStratificationDimension)

Maximum stratification dimension accessor.

Parameters
maxint

Maximum stratification dimension.

setName(name)

Accessor to the object’s name.

Parameters
namestr

The name of the object.

setPartialStratification(partialStratification)

Partial stratification accessor.

Parameters
partialStratificationbool

Partial stratification.

setProgressCallback(*args)

Set up a progress callback.

Can be used to programmatically report the progress of a simulation.

Parameters
callbackcallable

Takes a float as argument as percentage of progress.

Examples

>>> import sys
>>> import openturns as ot
>>> experiment = ot.MonteCarloExperiment()
>>> X = ot.RandomVector(ot.Normal())
>>> Y = ot.CompositeRandomVector(ot.SymbolicFunction(['X'], ['1.1*X']), X)
>>> event = ot.ThresholdEvent(Y, ot.Less(), -2.0)
>>> algo = ot.ProbabilitySimulationAlgorithm(event, experiment)
>>> algo.setMaximumOuterSampling(100)
>>> algo.setMaximumCoefficientOfVariation(-1.0)
>>> def report_progress(progress):
...     sys.stderr.write('-- progress=' + str(progress) + '%\n')
>>> algo.setProgressCallback(report_progress)
>>> algo.run()
setQuadrantOrientation(quadrantOrientation)

Quadrant orientation accessor.

Parameters
orientationsequence of float

Quadrant orientation.

setRootStrategy(rootStrategy)

Set the root strategy.

Parameters
strategyRootStrategy

Root strategy adopted.

setSamplingStrategy(samplingStrategy)

Set the direction sampling strategy.

Parameters
strategySamplingStrategy

Direction sampling strategy adopted.

setShadowedId(id)

Accessor to the object’s shadowed id.

Parameters
idint

Internal unique identifier.

setStopCallback(*args)

Set up a stop callback.

Can be used to programmatically stop a simulation.

Parameters
callbackcallable

Returns an int deciding whether to stop or continue.

Examples

Stop a Monte Carlo simulation algorithm using a time limit

>>> import openturns as ot
>>> experiment = ot.MonteCarloExperiment()
>>> X = ot.RandomVector(ot.Normal())
>>> Y = ot.CompositeRandomVector(ot.SymbolicFunction(['X'], ['1.1*X']), X)
>>> event = ot.ThresholdEvent(Y, ot.Less(), -2.0)
>>> algo = ot.ProbabilitySimulationAlgorithm(event, experiment)
>>> algo.setMaximumOuterSampling(10000000)
>>> algo.setMaximumCoefficientOfVariation(-1.0)
>>> timer = ot.TimerCallback(0.1)
>>> algo.setStopCallback(timer)
>>> algo.run()
setVerbose(verbose)

Accessor to verbosity.

Parameters
verbosity_enabledbool

If True, make the computation verbose. By default it is verbose.

setVisibility(visible)

Accessor to the object’s visibility state.

Parameters
visiblebool

Visibility flag.