First order second moment

ffpack.rrm.firstOrderSecondMoment.mvalFOSM(dim, g, dg, mus, sigmas, dx=1e-06)

First order second moment method based on mean value algorithm.

Parameters
  • dim (integer) – Space dimension ( number of random variables ).

  • g (function) – Limit state function. It will be called like g( [ x1, x2, … ] ).

  • dg (array_like of function) – Gradient of the limit state function. It should be an array_like of function like dg = [ dg_dx1, dg_dx2, … ]. To get the derivative of i-th random variable at ( x1*, x2*, … ), dg[ i ]( x1*, x2*, … ) will be called. dg can be None, see the following Notes.

  • mus (1d array) – Mean of the random variables.

  • sigmas (1d array) – Variance of the random variables.

  • dx (scalar, optional) – Spacing for auto differentiation. Not required if dg is provided.

Returns

  • beta (scalar) – Reliability index.

  • pf (scalar) – probability of failure.

Raises

ValueError – If the dim is less than 1. If the dim does not match the length of mus and sigmas.

Notes

If dg is None, the numerical differentiation will be used. The tolerance of the numerical differentiation can be changed in globalConfig.

Examples

>>> from ffpack.rrm import mvalFOSM
>>> dim = 2
>>> g = lambda X: 3 * X[ 0 ] - 2 * X[ 1 ]
>>> dg = [ lambda X: 3, lambda X: -2 ]
>>> mus = [ 1, 1 ]
>>> sigmas = [ 3, 4 ]
>>> beta, pf = mvalFOSM( dim, g, dg, mus, sigmas)

First order reliability method

ffpack.rrm.firstOrderReliabilityMethod.coptFORM(dim, g, distObjs, corrMat, quadDeg=99, quadRange=8)

First order reliability method based on constrained optimization.

Parameters
  • dim (integer) – Space dimension ( number of random variables ).

  • g (function) – Limit state function. It will be called like g( [ x1, x2, … ] ).

  • distObjs (array_like of distributions) – Marginal distribution objects. It should be the freezed distribution objects with pdf, cdf, ppf. We recommend to use scipy.stats functions.

  • corrMat (2d matrix) – Correlation matrix of the marginal distributions.

  • quadDeg (integer) – Quadrature degree for Nataf transformation

  • quadRange (scalar) – Quadrature range for Nataf transformation. The integral will be performed in the range [ -quadRange, quadRange ].

Returns

  • beta (scalar) – Reliability index.

  • pf (scalar) – Probability of failure.

  • uCoord (1d array) – Design point coordinate in U space.

  • xCoord (1d array) – Design point coordinate in X space.

Raises

ValueError – If the dim is less than 1. If the dim does not match the disObjs and corrMat. If corrMat is not 2d matrix. If corrMat is not positive definite. If corrMat is not symmetric. If corrMat diagonal is not 1.

Examples

>>> from ffpack.rrm import coptFORM
>>> dim = 2
>>> g = lambda X: -np.sum( X ) + 1
>>> distObjs = [ stats.norm(), stats.norm() ]
>>> corrMat = np.eye( dim )
>>> beta, pf, uCoord, xCoord = coptFORM( dim, g, distObjs, corrMat )
ffpack.rrm.firstOrderReliabilityMethod.hlrfFORM(dim, g, dg, distObjs, corrMat, iter=1000, tol=1e-06, quadDeg=99, quadRange=8, dx=1e-06)

First order reliability method based on Hasofer-Lind-Rackwitz-Fiessler algorithm.

Parameters
  • dim (integer) – Space dimension ( number of random variables ).

  • g (function) – Limit state function. It will be called like g( [ x1, x2, … ] ).

  • dg (array_like of function) – Gradient of the limit state function. It should be an array_like of function like dg = [ dg_dx1, dg_dx2, … ]. To get the derivative of i-th random variable at ( x1*, x2*, … ), dg[ i ]( x1*, x2*, … ) will be called. dg can be None, see the following Notes.

  • distObjs (array_like of distributions) – Marginal distribution objects. It should be the freezed distribution objects with pdf, cdf, ppf. We recommend to use scipy.stats functions.

  • corrMat (2d matrix) – Correlation matrix of the marginal distributions.

  • iter (integer) – Maximum iteration steps.

  • tol (scalar) – Tolerance to demtermine if the iteration converges.

  • quadDeg (integer) – Quadrature degree for Nataf transformation

  • quadRange (scalar) – Quadrature range for Nataf transformation. The integral will be performed in the range [ -quadRange, quadRange ].

  • dx (scalar, optional) – Spacing for auto differentiation. Not required if dg is provided.

Returns

  • beta (scalar) – Reliability index.

  • pf (scalar) – Probability of failure.

  • uCoord (1d array) – Design point coordinate in U space.

  • xCoord (1d array) – Design point coordinate in X space.

Raises

ValueError – If the dim is less than 1. If the dim does not match the disObjs and corrMat. If corrMat is not 2d matrix. If corrMat is not positive definite. If corrMat is not symmetric. If corrMat diagonal is not 1.

Notes

If dg is None, the numerical differentiation will be used. The tolerance of the numerical differentiation can be changed in globalConfig.

Examples

>>> from ffpack.rrm import hlrfFORM
>>> dim = 2
>>> g = lambda X: -np.sum( X ) + 1
>>> dg = [ lambda X: -1, lambda X: -1 ]
>>> distObjs = [ stats.norm(), stats.norm() ]
>>> corrMat = np.eye( dim )
>>> beta, pf, uCoord, xCoord = hlrfFORM( dim, g, dg, distObjs, corrMat )

Second order reliability method

ffpack.rrm.secondOrderReliabilityMethod.breitungSORM(dim, g, dg, distObjs, corrMat, quadDeg=99, quadRange=8, dx=1e-06)

Second order reliability method based on Breitung algorithm.

Parameters
  • dim (integer) – Space dimension ( number of random variables ).

  • g (function) – Limit state function. It will be called like g( [ x1, x2, … ] ).

  • dg (array_like of function) – Gradient of the limit state function. It should be an array_like of function like dg = [ dg_dx1, dg_dx2, … ]. To get the derivative of i-th random variable at ( x1*, x2*, … ), dg[ i ]( x1*, x2*, … ) will be called. dg can be None, see the following Notes.

  • distObjs (array_like of distributions) – Marginal distribution objects. It should be the freezed distribution objects with pdf, cdf, ppf. We recommend to use scipy.stats functions.

  • corrMat (2d matrix) – Correlation matrix of the marginal distributions.

  • quadDeg (integer) – Quadrature degree for Nataf transformation

  • quadRange (scalar) – Quadrature range for Nataf transformation. The integral will be performed in the range [ -quadRange, quadRange ].

  • dx (scalar, optional) – Spacing for auto differentiation. Not required if dg is provided.

Returns

  • beta (scalar) – Reliability index.

  • pf (scalar) – Probability of failure.

  • uCoord (1d array) – Design point coordinate in U space.

  • xCoord (1d array) – Design point coordinate in X space.

Raises

ValueError – If the dim is less than 1. If the dim does not match the disObjs and corrMat. If corrMat is not 2d matrix. If corrMat is not positive definite. If corrMat is not symmetric. If corrMat diagonal is not 1.

Notes

If dg is None, the numerical differentiation will be used. The tolerance of the numerical differentiation can be changed in globalConfig.

Examples

>>> from ffpack.rrm import breitungSORM
>>> dim = 2
>>> g = lambda X: -np.sum( X ) + 1
>>> dg = [ lambda X: -1, lambda X: -1 ]
>>> distObjs = [ stats.norm(), stats.norm() ]
>>> corrMat = np.eye( dim )
>>> beta, pf, uCoord, xCoord = breitungSORM( dim, g, dg, distObjs, corrMat )
ffpack.rrm.secondOrderReliabilityMethod.hrackSORM(dim, g, dg, distObjs, corrMat, quadDeg=99, quadRange=8, dx=1e-06)

Second order reliability method based on Hohenbichler and Rackwitz algorithm.

Parameters
  • dim (integer) – Space dimension ( number of random variables ).

  • g (function) – Limit state function. It will be called like g( [ x1, x2, … ] ).

  • dg (array_like of function) – Gradient of the limit state function. It should be an array_like of function like dg = [ dg_dx1, dg_dx2, … ]. To get the derivative of i-th random variable at ( x1*, x2*, … ), dg[ i ]( x1*, x2*, … ) will be called. dg can be None, see the following Notes.

  • distObjs (array_like of distributions) – Marginal distribution objects. It should be the freezed distribution objects with pdf, cdf, ppf. We recommend to use scipy.stats functions.

  • corrMat (2d matrix) – Correlation matrix of the marginal distributions.

  • quadDeg (integer) – Quadrature degree for Nataf transformation

  • quadRange (scalar) – Quadrature range for Nataf transformation. The integral will be performed in the range [ -quadRange, quadRange ].

  • dx (scalar, optional) – Spacing for auto differentiation. Not required if dg is provided.

Returns

  • beta (scalar) – Reliability index.

  • pf (scalar) – Probability of failure.

  • uCoord (1d array) – Design point coordinate in U space.

  • xCoord (1d array) – Design point coordinate in X space.

Raises

ValueError – If the dim is less than 1. If the dim does not match the disObjs and corrMat. If corrMat is not 2d matrix. If corrMat is not positive definite. If corrMat is not symmetric. If corrMat diagonal is not 1.

Notes

If dg is None, the numerical differentiation will be used. The tolerance of the numerical differentiation can be changed in globalConfig.

Examples

>>> from ffpack.rrm import tvedtSORM
>>> dim = 2
>>> g = lambda X: -np.sum( X ) + 1
>>> dg = [ lambda X: -1, lambda X: -1 ]
>>> distObjs = [ stats.norm(), stats.norm() ]
>>> corrMat = np.eye( dim )
>>> beta, pf, uCoord, xCoord = hrackSORM( dim, g, dg, distObjs, corrMat )
ffpack.rrm.secondOrderReliabilityMethod.tvedtSORM(dim, g, dg, distObjs, corrMat, quadDeg=99, quadRange=8, dx=1e-06)

Second order reliability method based on Tvedt algorithm.

Parameters
  • dim (integer) – Space dimension ( number of random variables ).

  • g (function) – Limit state function. It will be called like g( [ x1, x2, … ] ).

  • dg (array_like of function) – Gradient of the limit state function. It should be an array_like of function like dg = [ dg_dx1, dg_dx2, … ]. To get the derivative of i-th random variable at ( x1*, x2*, … ), dg[ i ]( x1*, x2*, … ) will be called. dg can be None, see the following Notes.

  • distObjs (array_like of distributions) – Marginal distribution objects. It should be the freezed distribution objects with pdf, cdf, ppf. We recommend to use scipy.stats functions.

  • corrMat (2d matrix) – Correlation matrix of the marginal distributions.

  • quadDeg (integer) – Quadrature degree for Nataf transformation

  • quadRange (scalar) – Quadrature range for Nataf transformation. The integral will be performed in the range [ -quadRange, quadRange ].

  • dx (scalar, optional) – Spacing for auto differentiation. Not required if dg is provided.

Returns

  • beta (scalar) – Reliability index.

  • pf (scalar) – Probability of failure.

  • uCoord (1d array) – Design point coordinate in U space.

  • xCoord (1d array) – Design point coordinate in X space.

Raises

ValueError – If the dim is less than 1. If the dim does not match the disObjs and corrMat. If corrMat is not 2d matrix. If corrMat is not positive definite. If corrMat is not symmetric. If corrMat diagonal is not 1.

Notes

If dg is None, the numerical differentiation will be used. The tolerance of the numerical differentiation can be changed in globalConfig.

Examples

>>> from ffpack.rrm import tvedtSORM
>>> dim = 2
>>> g = lambda X: -np.sum( X ) + 1
>>> dg = [ lambda X: -1, lambda X: -1 ]
>>> distObjs = [ stats.norm(), stats.norm() ]
>>> corrMat = np.eye( dim )
>>> beta, pf, uCoord, xCoord = tvedtSORM( dim, g, dg, distObjs, corrMat )

Simulation based reliability method

ffpack.rrm.simulationBasedReliabilityMethod.subsetSimulation(dim, g, distObjs, corrMat, numSamples, maxSubsets, probLevel=0.1, quadDeg=99, quadRange=8, randomSeed=None)

Second order reliability method based on Breitung algorithm.

Parameters
  • dim (integer) – Space dimension ( number of random variables ).

  • g (function) – Limit state function. It will be called like g( [ x1, x2, … ] ).

  • distObjs (array_like of distributions) – Marginal distribution objects. It should be the freezed distribution objects with pdf, cdf, ppf. We recommend to use scipy.stats functions.

  • corrMat (2d matrix) – Correlation matrix of the marginal distributions.

  • numSamples (integer) – Number of samples in each subset.

  • maxSubsets (scalar) – Maximum number of subsets used to compute the failure probability.

  • probLevel (scalar, optional) – Probability level for intermediate subsets.

  • quadDeg (integer, optional) – Quadrature degree for Nataf transformation

  • quadRange (scalar, optional) – Quadrature range for Nataf transformation. The integral will be performed in the range [ -quadRange, quadRange ].

  • randomSeed (integer, optional) – Random seed. If randomSeed is none or is not an integer, the random seed in global config will be used.

Returns

  • pf (scalar) – Probability of failure.

  • allLsfValue (array_like) – Values of limit state function in each subset.

  • allUSamples (array_like) – Samples of U space in each subset.

  • allXSamples (array_like) – Samples of X space in each subset.

Raises

ValueError – If the dim is less than 1. If the dim does not match the disObjs and corrMat. If corrMat is not 2d matrix. If corrMat is not positive definite. If corrMat is not symmetric. If corrMat diagonal is not 1.

Notes

Nataf transformation is used for the marginal distributions.

Examples

>>> from ffpack.rrm import subsetSimulation
>>> dim = 2
>>> g = lambda X: -np.sum( X ) + 1
>>> distObjs = [ stats.norm(), stats.norm() ]
>>> corrMat = np.eye( dim )
>>> numSamples, maxSubsets = 500, 10
>>> pf = subsetSimulation( dim, g, distObjs, corrMat, numSamples, maxSubsets )