bet.sensitivity package

Submodules

bet.sensitivity.chooseQoIs module

This module contains functions for choosing optimal sets of QoIs to use in the stochastic inverse problem.

bet.sensitivity.chooseQoIs.calculate_avg_condnum(input_set, qoi_set=None)

Given gradient vectors at some points (centers) in the input space and given a specific set of QoIs, caculate the average condition number of the matrices formed by the gradient vectors of each QoI map at each center.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None

  • qoi_set (list) – List of QoI indices

Return type

tuple

Returns

(condnum, singvals) where condnum is a float and singvals has shape (num_centers, output_dim)

bet.sensitivity.chooseQoIs.calculate_avg_measure(input_set, qoi_set=None, bin_measure=None)

If you are using bin_ratio to define the hyperrectangle in the output space you must give this method gradient vectors normalized with respect to the 1-norm. If you are using bin_size to define the hyperrectangle in the output space you must give this method the original gradient vectors. If you also give a bin_measure, this method will approximate the measure of the region of non-zero probability in the inverse solution. Given gradient vectors at some points (centers) in the input space and given a specific set of QoIs, calculate the expected measure of the inverse image of a box in the data space using local linear approximations of the map Q.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None

  • qoi_set (list) – List of QoI indices

  • bin_measure (float) – The measure of the output_dim hyperrectangle to invert into the input space

Return type

tuple

Returns

(avg_measure, singvals) where avg_measure is a float and singvals has shape (num_centers, output_dim)

bet.sensitivity.chooseQoIs.calculate_avg_skewness(input_set, qoi_set=None)

Given gradient vectors at some points (centers) in the input space and given a specific set of QoIs, caculate the average skewness of the arrays formed by the gradient vectors of each QoI map at each center.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None

  • qoi_set (list) – List of QoI indices

Return type

tuple

Returns

(hmean_skewG, skewgi) where hmean_skewG is the harmonic mean of skewness at each center in the input space (float) and skewgi has shape (num_centers, output_dim)

bet.sensitivity.chooseQoIs.chooseOptQoIs(input_set, qoiIndices=None, num_qois_return=None, num_optsets_return=None, inner_prod_tol=1.0, measure=False, remove_zeros=True)

Given gradient vectors at some points (centers) in the parameter space, a set of QoIs to choose from, and the number of desired QoIs to return, this method returns the num_optsets_return best sets of QoIs with with repsect to either the average measure of the matrix formed by the gradient vectors of each QoI map, OR the average skewness of the inverse image of this set of QoIs, computed as the product of the singular values of the same matrix. This method is brute force, i.e., if the method is given 10,000 QoIs and told to return the N best sets of 3, it will check all 10,000 choose 3 possible sets. See chooseOptQoIs_large for a less computationally expensive approach.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None

  • qoiIndices (np.ndarray of size (1, num QoIs to consider)) – Set of QoIs to consider. Default is xrange(0, input_set._jacobians.shape[1])

  • num_qois_return (int) – Number of desired QoIs to use in the inverse problem. Default is input_dim

  • num_optsets_return (int) – Number of best sets to return Default is 10

  • measure (boolean) – If measure is True, use calculate_avg_measure to determine optimal QoIs, else use calculate_avg_skewness

  • remove_zeros (boolean) – If True, find_unique_vecs will remove any QoIs that have a zero gradient

Return type

np.ndarray of shape (num_optsets_returned, num_qois_returned + 1)

Returns

measure_skewness_indices_mat

bet.sensitivity.chooseQoIs.chooseOptQoIs_large(input_set, qoiIndices=None, max_qois_return=None, num_optsets_return=None, inner_prod_tol=None, measskew_tol=None, measure=False, remove_zeros=True)

Given gradient vectors at some points (centers) in the input space, a large set of QoIs to choose from, and the number of desired QoIs to return, this method returns the set of optimal QoIs of size 2, 3, … max_qois_return to use in the inverse problem by choosing the sets with the smallest average measure(skewness).

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None.

  • qoiIndices (np.ndarray of size (1, num QoIs to consider)) – Set of QoIs to consider from input_set._jacobians. Default is xrange(0, input_set._jacobians.shape[1])

  • max_qois_return (int) – Maximum number of desired QoIs to use in the inverse problem. Default is input_dim

  • num_optsets_return (int) – Number of best sets to return Default is 10

  • inner_prod_tol (float) – Maximum acceptable average inner product between two QoI maps.

  • measskew_tol (float) – Throw out all sets of QoIs with average measure(skewness) number greater than this.

  • measure (boolean) – If measure is True, use calculate_avg_measure to determine optimal QoIs, else use calculate_avg_skewness

  • remove_zeros (boolean) – If True, find_unique_vecs will remove any QoIs that have a zero gradient vector at atleast one point in \(\Lambda\).

Return type

tuple

Returns

(measure_skewness_indices_mat, optsingvals) where measure_skewness_indices_mat has shape (num_optsets_return, num_qois_return+1) and optsingvals has shape (num_centers, num_qois_return, num_optsets_return)

bet.sensitivity.chooseQoIs.chooseOptQoIs_large_verbose(input_set, qoiIndices=None, max_qois_return=None, num_optsets_return=None, inner_prod_tol=None, measskew_tol=None, measure=False, remove_zeros=True)

Given gradient vectors at some points (centers) in the parameter space, a large set of QoIs to choose from, and the number of desired QoIs to return, this method return the set of optimal QoIs of size 1, 2, … max_qois_return to use in the inverse problem by choosing the set with smallext average skewness. Also a tensor that represents the singular values of the matrices formed by the gradient vectors of the optimal QoIs at each center is returned.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None.

  • qoiIndices (np.ndarray of size (1, num QoIs to consider)) – Set of QoIs to consider from G. Default is xrange(0, G.shape[1]).

  • max_qois_return (int) – Maximum number of desired QoIs to use in the inverse problem. Default is input_dim.

  • num_optsets_return (int) – Number of best sets to return. Default is 10.

  • inner_prod_tol (float) – Throw out one vectors from each pair of QoIs that has average inner product greater than this. Default is 0.9.

  • measskew_tol (float) – Throw out all sets of QoIs with average measure(skewness) number greater than this. Default is max_float.

  • measure (boolean) – If measure is True, use calculate_avg_measure to determine optimal QoIs, else use calculate_avg_skewness

  • remove_zeros (boolean) – If True, find_unique_vecs will remove any QoIs that have a zero gradient vector at atleast one point in \(\Lambda\).

Return type

tuple

Returns

(measure_skewness_indices_mat, optsingvals) where measure_skewness_indices_mat has shape (num_optsets_return, num_qois_return+1) and optsingvals is a list where each element has shape (num_centers, num_qois_return, num_optsets_return). num_qois_return will change for each element of the list.

bet.sensitivity.chooseQoIs.chooseOptQoIs_verbose(input_set, qoiIndices=None, num_qois_return=None, num_optsets_return=None, inner_prod_tol=1.0, measure=False, remove_zeros=True)

Given gradient vectors at some points (centers) in the parameter space, a set of QoIs to choose from, and the number of desired QoIs to return, this method returns the num_optsets_return best sets of QoIs with with repsect to either the average measure of the matrix formed by the gradient vectors of each QoI map, OR the average skewness of the inverse image of this set of QoIs, computed as the product of the singular values of the same matrix. This method is brute force, i.e., if the method is given 10,000 QoIs and told to return the N best sets of 3, it will check all 10,000 choose 3 possible sets. See chooseOptQoIs_large for a less computationally expensive approach.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None

  • qoiIndices (np.ndarray of size (1, num QoIs to consider)) – Set of QoIs to consider. Default is xrange(0, input_set._jacobians.shape[1])

  • num_qois_return (int) – Number of desired QoIs to use in the inverse problem. Default is input_dim

  • num_optsets_return (int) – Number of best sets to return Default is 10

  • measure (boolean) – If measure is True, use calculate_avg_measure to determine optimal QoIs, else use calculate_avg_skewness

  • remove_zeros (boolean) – If True, find_unique_vecs will remove any QoIs that have a zero gradient

Return type

np.ndarray of shape (num_optsets_returned, num_qois_returned + 1)

Returns

measure_skewness_indices_mat

bet.sensitivity.chooseQoIs.find_good_sets(input_set, good_sets_prev, unique_indices, num_optsets_return, measskew_tol, measure)

Todo

Use the idea we only know vectors are with 10% accuracy to guide inner_prod tol and skewness_tol.

Given gradient vectors at each center in the parameter space and given good sets of size (n - 1), return good sets of size n. That is, return sets of size n that have average measure(skewness) less than some tolerance.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None.

  • good_sets_prev (np.ndarray of size (num_good_sets_prev, n - 1)) – Good sets of QoIs of size n - 1.

  • unique_indices (np.ndarray of size (num_unique_qois, 1)) – Unique QoIs to consider.

  • num_optsets_return (int) – Number of best sets to return

  • measskew_tol (float) – Throw out all sets of QoIs with average measure(skewness) number greater than this.

  • measure (boolean) – If measure is True, use calculate_avg_measure to determine optimal QoIs, else use calculate_avg_skewness

Return type

tuple

Returns

(good_sets, best_sets, optsingvals_tensor) where good sets has size (num_good_sets, n), best sets has size (num_optsets_return, n + 1) and optsingvals_tensor has size (num_centers, n, input_dim)

bet.sensitivity.chooseQoIs.find_unique_vecs(input_set, inner_prod_tol, qoiIndices=None, remove_zeros=True)

Given gradient vectors at each center in the parameter space, sort throught them and remove any QoI that has a zero vector at any center, then remove one from any pair of QoIs that have an average inner product greater than some tolerance, i.e., an average angle between the two vectors smaller than some tolerance.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _jacobians is not None

  • inner_prod_tol (float) – Maximum acceptable average inner product between two QoI maps

  • qoiIndices (:class:’np.ndarray of size (1, num QoIs to consider)) – Set of QoIs to consider

  • remove_zeros (boolean) – If True, find_unique_vecs will remove any QoIs that have a zero gradient vector at atleast one point in \(\Lambda\)

Return type

np.ndarray of shape (num_unique_vecs, 1)

Returns

unique_vecs

bet.sensitivity.gradients module

This module contains functions for approximating jacobians of QoI maps. All methods that cluster points around centers are written to return the input_set._values in the following order : CENTERS, FOLLOWED BY THE CLUSTER AROUND THE FIRST CENTER, THEN THE CLUSTER AROUND THE SECOND CENTER AND SO ON.

bet.sensitivity.gradients.calculate_gradients_cfd(cluster_discretization, normalize=True)

Approximate gradient vectors at num_centers, centers.shape[0] points in the parameter space for each QoI map. THIS METHOD IS DEPENDENT ON USING :meth:~bet.sensitivity.pick_cfd_points TO CHOOSE SAMPLES FOR THE CFD STENCIL AROUND EACH CENTER. THE ORDERING MATTERS.

Parameters
  • cluster_discretization (discretization) – Must contain input and output values for the sample clusters.

  • normalize (boolean) – If normalize is True, normalize each gradient vector

Return type

discretization

Returns

A new discretization that contains only the centers of the clusters and their associated _jacobians which are tensor representation of the gradient vectors of each QoI map at each point in centers numpy.ndarray of shape (num_samples, output_dim, input_dim)

bet.sensitivity.gradients.calculate_gradients_ffd(cluster_discretization, normalize=True)

Approximate gradient vectors at num_centers, centers.shape[0] points in the parameter space for each QoI map. THIS METHOD IS DEPENDENT ON USING :meth:~bet.sensitivity.gradients.pick_ffd_points TO CHOOSE SAMPLES FOR THE FFD STENCIL AROUND EACH CENTER. THE ORDERING MATTERS.

Parameters
  • cluster_discretization (discretization) – Must contain input and output values for the sample clusters.

  • normalize (boolean) – If normalize is True, normalize each gradient vector

Return type

discretization

Returns

A new discretization that contains only the centers of the clusters and their associated _jacobians which are tensor representation of the gradient vectors of each QoI map at each point in centers numpy.ndarray of shape (num_samples, output_dim, input_dim)

bet.sensitivity.gradients.calculate_gradients_rbf(cluster_discretization, num_centers=None, num_neighbors=None, RBF=None, ep=None, normalize=True)

Approximate gradient vectors at num_centers, centers.shape[0] points in the parameter space for each QoI map using a radial basis function interpolation method.

Parameters
  • cluster_discretization (discretization) – Must contain input and output values for the sample clusters.

  • num_centers (int) – The number of cluster centers.

  • num_neighbors (int) – Number of nearest neighbors to use in gradient approximation. Default value is input_dim + 2

  • RBF (string) – Choice of radial basis function. Default is Gaussian

  • ep (float) – Choice of shape parameter for radial basis function. Default value is 1.0

  • normalize (boolean) – If normalize is True, normalize each gradient vector

Return type

discretization

Returns

A new discretization that contains only the centers of the clusters and their associated _jacobians which are tensor representation of the gradient vectors of each QoI map at each point in centers numpy.ndarray of shape (num_samples, output_dim, input_dim)

bet.sensitivity.gradients.pick_cfd_points(input_set, radii_vec)

Pick 2*input_dim points, for each center, for centered finite difference gradient approximation. The center are not needed for the CFD gradient approximation, they are returned for consistency with the other methods and because of the common need to have not just the gradient but also the QoI value at the centers in adaptive sampling algorithms.The points are returned in the order: centers, followed by the cluster around the first center, then the cluster around the second center and so on.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _values is not None

  • radii_vec (numpy.ndarray of shape (input_dim,)) – The radius of the stencil, along each axis

Return type

sample_set

Returns

Centers and clusters of samples near each center (values are numpy.ndarray of shape ((num_close+1)*``num_centers``, input_dim))

bet.sensitivity.gradients.pick_ffd_points(input_set, radii_vec)

Pick input_dim points, for each centers, for a forward finite difference gradient approximation. The points are returned in the order: centers, followed by the cluster around the first center, then the cluster around the second center and so on.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _values is not None

  • radii_vec (numpy.ndarray of shape (input_dim,)) – The radius of the stencil, along each axis

Return type

sample_set

Returns

Centers and clusters of samples near each center (values are numpy.ndarray of shape ((num_close+1)*``num_centers``, input_dim))

bet.sensitivity.gradients.radial_basis_function(r, kernel=None, ep=None)

Evaluate a chosen radial basis function. Allow for the choice of several radial basis functions to use in :meth:~bet.sensitivity.gradients.calculate_gradients_rbf

Parameters
  • r (numpy.ndarray) – Distances from the reference point

  • kernel (string) – Choice of radial basis funtion. Default is C4Matern

  • ep (float) – Shape parameter for the radial basis function. Default is 1.0

Return type

numpy.ndarray of shape (r.shape)

Returns

Radial basis function evaluated for each element of r

bet.sensitivity.gradients.radial_basis_function_dxi(r, xi, kernel=None, ep=None)

Evaluate a partial derivative of a chosen radial basis function. Allow for the choice of several radial basis functions to use in the :meth:~bet.sensitivity.gradients.calculate_gradients_rbf.

Parameters
  • r (numpy.ndarray) – Distances from the reference point

  • xi (numpy.ndarray) – Distances from the reference point in dimension i

  • kernel (string) – Choice of radial basis funtion. Default is C4Matern

  • ep (float) – Shape parameter for the radial basis function. Default is 1.0

Return type

numpy.ndarray of shape (r.shape)

Returns

Radial basis function evaluated for each element of r

bet.sensitivity.gradients.sample_l1_ball(input_set, num_close, radii_vec)

Pick num_close points in a the L-1 ball of length 2*``radii_vec`` around a point in the input space, do this for each point in centers. If this box extends outside of the domain of the input space, we sample the intersection.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _values is not None

  • num_close (int) – Number of points in each cluster

  • radii_vec (numpy.ndarray of shape (input_dim,)) – Each side of the box will have length 2*radii_vec[i]

Return type

sample_set

Returns

Centers and clusters of samples near each center (values are numpy.ndarray of shape ((num_close+1)*``num_centers``, input_dim))

bet.sensitivity.gradients.sample_linf_ball(input_set, num_close, radii_vec)

Pick num_close points in a the L-inifity ball of length 2*``radii_vec`` around a point in the input space, do this for each point in centers. If this box extends outside of the domain of the input space, we sample the intersection.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _values is not None

  • num_close (int) – Number of points in each cluster

  • radii_vec (numpy.ndarray of shape (input_dim,)) – Each side of the box will have length 2*radii_vec[i]

Return type

sample_set

Returns

Centers and clusters of samples near each center (values are numpy.ndarray of shape ((num_close+1)*``num_centers``, input_dim))

bet.sensitivity.gradients.sample_lp_ball(input_set, num_close, radius, p_num=2)

Pick num_close points in a the Lp ball of length 2*``radii_vec`` around a point in the input space, do this for each point in centers. If this box extends outside of the domain of the input space, we sample the intersection.

Parameters
  • input_set (sample_set) – The input sample set. Make sure the attribute _values is not None

  • num_close (int) – Number of points in each cluster

  • radii_vec (numpy.ndarray of shape (input_dim,)) – Each side of the box will have length 2*radii_vec[i]

  • p_num (float) – \(0 < p \leq \infty\), p for the lp norm where infinity is numpy.inf

Return type

sample_set

Returns

Centers and clusters of samples near each center (values are numpy.ndarray of shape ((num_close+1)*``num_centers``, input_dim))

Module contents

This subpackage provides methods for approximating gradients of QoI maps and choosing optimal QoIs to use in the inverse problem.

  • gradients provides methods for approximating gradients of QoI maps.

  • chooseQoIs provides methods for choosing optimal QoIs to use in the inverse problem.