ThresholdEvent

class ThresholdEvent(*args)

Random vector defined from a comparison operator and a threshold.

The event occurs when the realization of the underlying random vector exceeds the threshold.

Parameters:
antecedentRandomVector of dimension 1

Output variable of interest.

comparisonOperatorComparisonOperator

Comparison operator used to compare antecedent with threshold.

thresholdfloat

threshold we want to compare to antecedent.

Notes

An event is defined as follows:

\cD_f = \{\vect{X} \in \Rset^n \, / \, g(\vect{X},\vect{d}) \le 0\}

where \vect{X} denotes a random input vector, representing the sources of uncertainties, \vect{d} is a determinist vector, representing the fixed variables and g(\vect{X},\vect{d}) is the limit state function of the model. The probability content of the event \cD_f is P_f:

P_f = \int_{g(\vect{X},\vect{d})\le 0}f_\vect{X}(\vect{x})\di{\vect{x}}

Here, the event considered is explicited directly from the limit state function g(\vect{X}\,,\,\vect{d}) : this is the classical structural reliability formulation. However, if the event is a threshold exceedance, it is useful to explicit the variable of interest Z=\tilde{g}(\vect{X}\,,\,\vect{d}), evaluated from the model \tilde{g}(.). In that case, the event considered, associated to the threshold z_s has the formulation:

\cD_f = \{ \vect{X} \in \Rset^n \, / \, Z=\tilde{g}(\vect{X}\,,\,\vect{d}) > z_s \}

and the limit state function is:

g(\vect{X}\,,\,\vect{d}) &= z_s - Z \\
                         &= z_s - \tilde{g}(\vect{X}\,,\,\vect{d})

P_f is the threshold exceedance probability, defined as:

P_f &= P(Z \geq z_s) \\
    &= \int_{g(\vect{X}, \vect{d}) \le 0} \pdf\di{\vect{x}}

Examples

>>> import openturns as ot
>>> dim = 2
>>> X = ot.RandomVector(ot.Normal(dim))
>>> model = ot.SymbolicFunction(['x1', 'x2'], ['x1+x2'])
>>> Y = ot.CompositeRandomVector(model, X)
>>> event = ot.ThresholdEvent(Y, ot.Less(), 1.0)

Methods

getAntecedent()

Accessor to the antecedent RandomVector in case of a composite RandomVector.

getClassName()

Accessor to the object's name.

getCovariance()

Accessor to the covariance of the RandomVector.

getDescription()

Accessor to the description of the RandomVector.

getDimension()

Accessor to the dimension of the RandomVector.

getDistribution()

Accessor to the distribution of the RandomVector.

getDomain()

Accessor to the domain of the Event.

getFrozenRealization(fixedPoint)

Compute realizations of the RandomVector.

getFrozenSample(fixedSample)

Compute realizations of the RandomVector.

getFunction()

Accessor to the Function in case of a composite RandomVector.

getId()

Accessor to the object's id.

getImplementation()

Accessor to the underlying implementation.

getMarginal(*args)

Get the random vector corresponding to the i^{th} marginal component(s).

getMean()

Accessor to the mean of the RandomVector.

getName()

Accessor to the object's name.

getOperator()

Accessor to the comparaison operator of the Event.

getParameter()

Accessor to the parameter of the distribution.

getParameterDescription()

Accessor to the parameter description of the distribution.

getRealization()

Compute one realization of the RandomVector.

getSample(size)

Compute realizations of the RandomVector.

getThreshold()

Accessor to the threshold of the Event.

intersect(other)

Intersection of two events.

isComposite()

Accessor to know if the RandomVector is a composite one.

isEvent()

Whether the random vector is an event.

join(other)

Union of two events.

setDescription(description)

Accessor to the description of the RandomVector.

setName(name)

Accessor to the object's name.

setParameter(parameters)

Accessor to the parameter of the distribution.

__init__(*args)
getAntecedent()

Accessor to the antecedent RandomVector in case of a composite RandomVector.

Returns:
antecedentRandomVector

Antecedent RandomVector \vect{X} in case of a CompositeRandomVector such as: \vect{Y}=f(\vect{X}).

getClassName()

Accessor to the object’s name.

Returns:
class_namestr

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

getCovariance()

Accessor to the covariance of the RandomVector.

Returns:
covarianceCovarianceMatrix

Covariance of the considered UsualRandomVector.

Examples

>>> import openturns as ot
>>> distribution = ot.Normal([0.0, 0.5], [1.0, 1.5], ot.CorrelationMatrix(2))
>>> randomVector = ot.RandomVector(distribution)
>>> ot.RandomGenerator.SetSeed(0)
>>> print(randomVector.getCovariance())
[[ 1    0    ]
 [ 0    2.25 ]]
getDescription()

Accessor to the description of the RandomVector.

Returns:
descriptionDescription

Describes the components of the RandomVector.

getDimension()

Accessor to the dimension of the RandomVector.

Returns:
dimensionpositive int

Dimension of the RandomVector.

getDistribution()

Accessor to the distribution of the RandomVector.

Returns:
distributionDistribution

Distribution of the considered UsualRandomVector.

Examples

>>> import openturns as ot
>>> distribution = ot.Normal([0.0, 0.0], [1.0, 1.0], ot.CorrelationMatrix(2))
>>> randomVector = ot.RandomVector(distribution)
>>> ot.RandomGenerator.SetSeed(0)
>>> print(randomVector.getDistribution())
Normal(mu = [0,0], sigma = [1,1], R = [[ 1 0 ]
 [ 0 1 ]])
getDomain()

Accessor to the domain of the Event.

Returns:
domainDomain

Describes the domain of an event.

getFrozenRealization(fixedPoint)

Compute realizations of the RandomVector.

In the case of a CompositeRandomVector or an event of some kind, this method returns the value taken by the random vector if the root cause takes the value given as argument.

Parameters:
fixedPointPoint

Point chosen as the root cause of the random vector.

Returns:
realizationPoint

The realization corresponding to the chosen root cause.

Examples

>>> import openturns as ot
>>> distribution = ot.Normal()
>>> randomVector = ot.RandomVector(distribution)
>>> f = ot.SymbolicFunction('x', 'x')
>>> compositeRandomVector = ot.CompositeRandomVector(f, randomVector)
>>> event = ot.ThresholdEvent(compositeRandomVector, ot.Less(), 0.0)
>>> print(event.getFrozenRealization([0.2]))
[0]
>>> print(event.getFrozenRealization([-0.1]))
[1]
getFrozenSample(fixedSample)

Compute realizations of the RandomVector.

In the case of a CompositeRandomVector or an event of some kind, this method returns the different values taken by the random vector when the root cause takes the values given as argument.

Parameters:
fixedSampleSample

Sample of root causes of the random vector.

Returns:
sampleSample

Sample of the realizations corresponding to the chosen root causes.

Examples

>>> import openturns as ot
>>> distribution = ot.Normal()
>>> randomVector = ot.RandomVector(distribution)
>>> f = ot.SymbolicFunction('x', 'x')
>>> compositeRandomVector = ot.CompositeRandomVector(f, randomVector)
>>> event = ot.ThresholdEvent(compositeRandomVector, ot.Less(), 0.0)
>>> print(event.getFrozenSample([[0.2], [-0.1]]))
    [ y0 ]
0 : [ 0  ]
1 : [ 1  ]
getFunction()

Accessor to the Function in case of a composite RandomVector.

Returns:
functionFunction

Function used to define a CompositeRandomVector as the image through this function of the antecedent \vect{X}: \vect{Y}=f(\vect{X}).

getId()

Accessor to the object’s id.

Returns:
idint

Internal unique identifier.

getImplementation()

Accessor to the underlying implementation.

Returns:
implImplementation

A copy of the underlying implementation object.

getMarginal(*args)

Get the random vector corresponding to the i^{th} marginal component(s).

Parameters:
iint or list of ints, 0\leq i < dim

Indicates the component(s) concerned. dim is the dimension of the RandomVector.

Returns:
vectorRandomVector

RandomVector restricted to the concerned components.

Notes

Let’s note \vect{Y}=\Tr{(Y_1,\dots,Y_n)} a random vector and I \in [1,n] a set of indices. If \vect{Y} is a UsualRandomVector, the subvector is defined by \tilde{\vect{Y}}=\Tr{(Y_i)}_{i \in I}. If \vect{Y} is a CompositeRandomVector, defined by \vect{Y}=f(\vect{X}) with f=(f_1,\dots,f_n), f_i some scalar functions, the subvector is \tilde{\vect{Y}}=(f_i(\vect{X}))_{i \in I}.

Examples

>>> import openturns as ot
>>> distribution = ot.Normal([0.0, 0.0], [1.0, 1.0], ot.CorrelationMatrix(2))
>>> randomVector = ot.RandomVector(distribution)
>>> ot.RandomGenerator.SetSeed(0)
>>> print(randomVector.getMarginal(1).getRealization())
[0.608202]
>>> print(randomVector.getMarginal(1).getDistribution())
Normal(mu = 0, sigma = 1)
getMean()

Accessor to the mean of the RandomVector.

Returns:
meanPoint

Mean of the considered UsualRandomVector.

Examples

>>> import openturns as ot
>>> distribution = ot.Normal([0.0, 0.5], [1.0, 1.5], ot.CorrelationMatrix(2))
>>> randomVector = ot.RandomVector(distribution)
>>> ot.RandomGenerator.SetSeed(0)
>>> print(randomVector.getMean())
[0,0.5]
getName()

Accessor to the object’s name.

Returns:
namestr

The name of the object.

getOperator()

Accessor to the comparaison operator of the Event.

Returns:
operatorComparisonOperator

Comparaison operator used to define the RandomVector.

getParameter()

Accessor to the parameter of the distribution.

Returns:
parameterPoint

Parameter values.

getParameterDescription()

Accessor to the parameter description of the distribution.

Returns:
descriptionDescription

Parameter names.

getRealization()

Compute one realization of the RandomVector.

Returns:
realizationPoint

Sequence of values randomly determined from the RandomVector definition. In the case of an event: one realization of the event (considered as a Bernoulli variable) which is a boolean value (1 for the realization of the event and 0 else).

Examples

>>> import openturns as ot
>>> distribution = ot.Normal([0.0, 0.0], [1.0, 1.0], ot.CorrelationMatrix(2))
>>> randomVector = ot.RandomVector(distribution)
>>> ot.RandomGenerator.SetSeed(0)
>>> print(randomVector.getRealization())
[0.608202,-1.26617]
>>> print(randomVector.getRealization())
[-0.438266,1.20548]
getSample(size)

Compute realizations of the RandomVector.

Parameters:
nint, n \geq 0

Number of realizations needed.

Returns:
realizationsSample

n sequences of values randomly determined from the RandomVector definition. In the case of an event: n realizations of the event (considered as a Bernoulli variable) which are boolean values (1 for the realization of the event and 0 else).

Examples

>>> import openturns as ot
>>> distribution = ot.Normal([0.0, 0.0], [1.0, 1.0], ot.CorrelationMatrix(2))
>>> randomVector = ot.RandomVector(distribution)
>>> ot.RandomGenerator.SetSeed(0)
>>> print(randomVector.getSample(3))
    [ X0        X1        ]
0 : [  0.608202 -1.26617  ]
1 : [ -0.438266  1.20548  ]
2 : [ -2.18139   0.350042 ]
getThreshold()

Accessor to the threshold of the Event.

Returns:
thresholdfloat

Threshold of the RandomVector.

intersect(other)

Intersection of two events.

Parameters:
eventRandomVector

A composite event

Returns:
eventRandomVector

Intersection event

isComposite()

Accessor to know if the RandomVector is a composite one.

Returns:
isCompositebool

Indicates if the RandomVector is of type Composite or not.

isEvent()

Whether the random vector is an event.

Returns:
isEventbool

Whether it takes it values in {0, 1}.

join(other)

Union of two events.

Parameters:
eventRandomVector

A composite event

Returns:
eventRandomVector

Union event

setDescription(description)

Accessor to the description of the RandomVector.

Parameters:
descriptionstr or sequence of str

Describes the components of the RandomVector.

setName(name)

Accessor to the object’s name.

Parameters:
namestr

The name of the object.

setParameter(parameters)

Accessor to the parameter of the distribution.

Parameters:
parametersequence of float

Parameter values.

Examples using the class

Simulate an Event

Simulate an Event

Estimate a probability with Monte Carlo

Estimate a probability with Monte Carlo

Use a randomized QMC algorithm

Use a randomized QMC algorithm

Use the Adaptive Directional Stratification Algorithm

Use the Adaptive Directional Stratification Algorithm

Use the post-analytical importance sampling algorithm

Use the post-analytical importance sampling algorithm

Use the Directional Sampling Algorithm

Use the Directional Sampling Algorithm

Create a threshold event

Create a threshold event

Specify a simulation algorithm

Specify a simulation algorithm

Estimate a flooding probability

Estimate a flooding probability

Use the Importance Sampling algorithm

Use the Importance Sampling algorithm

Estimate a probability with Monte-Carlo on axial stressed beam: a quick start guide to reliability

Estimate a probability with Monte-Carlo on axial stressed beam: a quick start guide to reliability

Estimate a buckling probability

Estimate a buckling probability

Exploitation of simulation algorithm results

Exploitation of simulation algorithm results

Use the FORM algorithm in case of several design points

Use the FORM algorithm in case of several design points

Subset Sampling

Subset Sampling

Use the FORM - SORM algorithms

Use the FORM - SORM algorithms

Non parametric Adaptive Importance Sampling (NAIS)

Non parametric Adaptive Importance Sampling (NAIS)

Test the design point with the Strong Maximum Test

Test the design point with the Strong Maximum Test

Time variant system reliability problem

Time variant system reliability problem

Create unions or intersections of events

Create unions or intersections of events

Axial stressed beam : comparing different methods to estimate a probability

Axial stressed beam : comparing different methods to estimate a probability

An illustrated example of a FORM probability estimate

An illustrated example of a FORM probability estimate

Cross Entropy Importance Sampling

Cross Entropy Importance Sampling

Using the FORM - SORM algorithms on a nonlinear function

Using the FORM - SORM algorithms on a nonlinear function

Control algorithm termination

Control algorithm termination