Active Subspaces

Active Subspaces module.

References
  • Paul Constantine. Active subspaces: Emerging ideas for dimension reduction in parameter studies, vol. 2 of SIAM Spotlights, SIAM, 2015.

  • Constantine et al. Python Active-subspaces Utility Library, Journal of Open Source Software, 1(5), 79, 2016.

ActiveSubspaces._bootstrap_replicate

Return a bootstrap replicate.

ActiveSubspaces._build_decompose_cov_matrix

Build and decompose the covariance matrix of the gradients.

ActiveSubspaces._compute_A_b

Compute the matrix A and the vector b to build a box around the inactive subspace for uniform sampling.

ActiveSubspaces._compute_bootstrap_ranges

Compute bootstrap ranges for eigenvalues and subspaces.

ActiveSubspaces._hit_and_run_inactive

A hit and run method for sampling the inactive variables from a polytope.

ActiveSubspaces._partition

Partition the eigenvectors to define the active and inactive subspaces.

ActiveSubspaces._rejection_sampling_inactive

A rejection sampling method for sampling the from a polytope.

ActiveSubspaces._rotate_x

A convenience function for rotating subspace coordinates to x space.

ActiveSubspaces._sample_inactive

Sample inactive variables.

ActiveSubspaces.activity_scores

Return the activity scores as defined in Constantine and Diaz, Global sensitivity metrics from active subspaces, arxiv.org/abs/1510.04361 Equation (21).

ActiveSubspaces.fit

Compute the active subspaces given the gradients of the model function wrt the input parameters, or given the input/outputs couples.

ActiveSubspaces.inverse_transform

Map the points in the active variable space to the original parameter space.

ActiveSubspaces.plot_eigenvalues

Plot the eigenvalues.

ActiveSubspaces.plot_eigenvectors

Plot the eigenvectors.

ActiveSubspaces.plot_sufficient_summary

Plot the sufficient summary.

ActiveSubspaces.transform

Map full variables to active and inactive variables.

class ActiveSubspaces(dim, method='exact', n_boot=100)[source]

Bases: athena.subspaces.Subspaces

Active Subspaces class

Parameters
  • dim (int) – dimension of the active subspace.

  • method (str) – method to compute the AS. Possible choices are ‘exact’ when the gradients are provided, or ‘local’ to use local linear models. This approach is related to the sufficient dimension reduction method known sometimes as the outer product of gradient method. See the 2001 paper ‘Structure adaptive approach for dimension reduction’ from Hristache, et al.

  • n_boot (int) – number of bootstrap samples. Default is 100.

_compute_A_b(reduced_input)[source]

Compute the matrix A and the vector b to build a box around the inactive subspace for uniform sampling.

Parameters

reduced_input (numpy.ndarray) – the value of the active variables.

Returns

matrix A, and vector b.

Return type

numpy.ndarray, numpy.ndarray

_frequent_directions(gradients)[source]

Function that performs the frequent directions algorithm for matrix sketching. For more details about the method, see “Frequent directions: Simple and deterministic matrix sketching.” Ghashami, Mina, et al. SIAM Journal on Computing 45.5 (2016): 1762-1792. doi: https://doi.org/10.1137/15M1009718

Parameters

gradients (iterable) – generator for spatial gradients.

Returns

the sorted eigenvalues and the corresponding eigenvectors for the reduced matrix.

Return type

numpy.ndarray, numpy.ndarray

_hit_and_run_inactive(reduced_input, n_points)[source]

A hit and run method for sampling the inactive variables from a polytope.

Parameters
  • reduced_input (numpy.ndarray) – the value of the active variables.

  • n_points (int) – the number of inactive variable samples.

Returns

n_points-by-(inactive_dim) matrix that contains values of the inactive variable that correspond to the given reduced_input.

Return type

numpy.ndarray

_rejection_sampling_inactive(reduced_input, n_points)[source]

A rejection sampling method for sampling the from a polytope.

Parameters
  • reduced_input (numpy.ndarray) – the value of the active variables.

  • n_points (int) – the number of inactive variable samples.

Returns

n_points-by-inactive_dim matrix that contains values of the inactive variable that correspond to the given reduced_input.

Return type

numpy.ndarray

_rotate_x(reduced_inputs, inactive_inputs)[source]

A convenience function for rotating subspace coordinates to x space.

Parameters
  • reduced_input (numpy.ndarray) – the value of the active variables.

  • inactive_inputs (numpy.ndarray) – the value of the inactive variables.

Returns

(n_samples * n_points)-by-n_params matrix that contains points in the original parameter space, (n_samples * n_points)-by-n_params matrix that contains integer indices. These indices identify which rows of the previous matrix (the full parameters) map to which rows of the active variables matrix.

Return type

numpy.ndarray, numpy.ndarray

_sample_inactive(reduced_input, n_points)[source]

Sample inactive variables.

Sample values of the inactive variables for a fixed value of the active variables when the original variables are bounded by a hypercube.

Parameters
  • reduced_input (numpy.ndarray) – the value of the active variables.

  • n_points (int) – the number of inactive variable samples,

Returns

n_points-by-(inactive_dim) matrix that contains values of the inactive variable that correspond to the given reduced_input.

Return type

numpy.ndarray

Note

The trick here is to sample the inactive variables z so that -1 <= W1*y + W2*z <= 1, where y is the given value of the active variables. In other words, we need to sample z such that it respects the linear inequalities W2*z <= 1 - W1*y, -W2*z <= 1 + W1*y. These inequalities define a polytope in R^(inactive_dim). We want to sample N points uniformly from the polytope. This function first tries a simple rejection sampling scheme, which finds a bounding hyperbox for the polytope, draws points uniformly from the bounding hyperbox, and rejects points outside the polytope. If that method does not return enough samples, the method tries a “hit and run” method for sampling from the polytope. If that does not work, it returns an array with N copies of a feasible point computed as the Chebyshev center of the polytope.

property activity_scores

Return the activity scores as defined in Constantine and Diaz, Global sensitivity metrics from active subspaces, arxiv.org/abs/1510.04361 Equation (21).

Returns

array with the activity score of each parameter.

Return type

numpy.ndarray

Raises

TypeError

Warning

self.fit has to be called in advance.

fit(inputs=None, outputs=None, gradients=None, weights=None, metric=None)[source]

Compute the active subspaces given the gradients of the model function wrt the input parameters, or given the input/outputs couples. Only two methods are available: ‘exact’ and ‘local’.

Parameters
  • inputs (numpy.ndarray) – input parameters oriented as rows.

  • outputs (numpy.ndarray) – corresponding outputs oriented as rows.

  • gradients (numpy.ndarray) – n_samples-by-n_params matrix containing the gradient samples oriented as rows. If frequent directions needed to be performed, gradients is an object of GeneratorType.

  • weights (numpy.ndarray) – n_samples-by-1 weight vector, corresponds to numerical quadrature rule used to estimate matrix whose eigenspaces define the active subspace.

  • metric (numpy.ndarray) – metric matrix output_dim-by-output-dim for vectorial active subspaces.

Raises

TypeError

Example
>>> # inputs shape is n_samples-by-n_params
>>> # outputs shape is n_samples-by-1
>>> # gradients shape is n_samples-by-n_params
>>> # if gradients are available use the 'exact' method:
>>> ss1 =  ActiveSubspaces(dim=1, method='exact', n_boot=150)
>>> ss1.fit(gradients=gradients)
>>> # for the frequent direction method to compute the eigenpairs
>>> # you need to pass a generator to gradients:
>>> gradients_gen = (grad for grad in gradients)
>>> ss2 =  ActiveSubspaces(dim=2, method='exact')
>>> ss2.fit(gradients=gradients_gen)
>>> # if no gradients are available use the 'local' method:
>>> ss3 =  ActiveSubspaces(dim=1, method='local', n_boot=150)
>>> ss3.fit(inputs=inputs, outputs=outputs)
inverse_transform(reduced_inputs, n_points=1)[source]

Map the points in the active variable space to the original parameter space.

Parameters
  • reduced_inputs (numpy.ndarray) – n_samples-by-dim matrix that contains points in the space of active variables.

  • n_points (int) – the number of points in the original parameter space that are returned that map to the given active variables. Defaults to 1.

Returns

(n_samples * n_points)-by-n_params matrix that contains points in the original parameter space, (n_samples * n_points)-by-n_params matrix that contains integer indices. These indices identify which rows of the previous matrix (the full parameters) map to which rows of the active variables matrix.

Return type

numpy.ndarray, numpy.ndarray

Raises

TypeError

Note

The inverse map depends critically on the self._sample_inactive method.

transform(inputs)[source]

Map full variables to active and inactive variables.

Points in the original input space are mapped to the active and inactive subspace.

Parameters

inputs (numpy.ndarray) – array n_samples-by-n_params containing the points in the original parameter space.

Returns

array n_samples-by-active_dim containing the mapped active variables; array n_samples-by-inactive_dim containing the mapped inactive variables.

Return type

numpy.ndarray, numpy.ndarray

Raises

TypeError