Feature Map¶
Module for the feature map class
- References
Francesco Romor, Marco Tezzele, Andrea Lario, Gianluigi Rozza. Kernel-based active subspaces with application to computational fluid dynamics parametric problems using discontinuous Galerkin method. International Journal for Numerical Methods in Engineering, 123(23):6000–6027, 2022. doi:10.1002/nme.7099
Sample the projection matrixx from the spectral distribution. |
|
Evaluate the feature map at inputs. |
|
Evaluate the Jacobian matrix of feature map at inputs. |
|
Get the projection matrix. |
|
Tune the parameters of the spectral distribution. |
|
Random Fourier Features map. |
|
Random Fourier Features map’s Jacobian. |
-
class
FeatureMap
(distr, bias, input_dim, n_features, params, sigma_f)[source] Bases:
object
Feature map class.
- Parameters
distr (str) – name of the spectral distribution to sample from the matrix.
bias (numpy.ndarray) – n_features dimensional bias. It is sampled from a Unifrom distribution on the interval [0, 2*PI].
input_dim (int) – dimension of the input space.
n_features (int) – dimension of the Reproducing Kernel Hilbert Space.
params (list) – number of hyperparameters of the spectral distribution.
sigma_f (int) – multiplicative constant of the feature map. Usually it is set as the empirical variance of the outputs.
- Variables
fmap (callable) – feature map used to project the input samples to the RKHS. Default value is rff_map.
fjac (callable) – Jacobian matrix of fmap. Default value is rff_jac, the analytical Jacobian of fmap.
_pr_matrix (numpy.ndarray) – n_features-by-input_dim projection matrix, whose rows are sampled from the spectral distribution distr.
:raises TypeError
-
_compute_pr_matrix
()[source] Sample the projection matrixx from the spectral distribution.
- Returns
the projection matrix.
- Return type
-
compute_fmap
(inputs)[source] Evaluate the feature map at inputs.
- Parameters
inputs (numpy.ndarray) – the inputs to project on the RKHS.
- Returns
the n_features dimensional projection of the inputs.
- Return type
-
compute_fmap_jac
(inputs)[source] Evaluate the Jacobian matrix of feature map at inputs.
- Parameters
inputs (numpy.ndarray) – the inputs at which compute the Jacobian matrix of the feature map.
- Returns
the n_features-by-input_dim dimensional Jacobian of the feature map at the inputs.
- Return type
-
property
pr_matrix
Get the projection matrix.
- Returns
the projection matrix.
- Return type
-
tune_pr_matrix
(func, bounds, method=None, maxiter=50, save_file=False, fn_args=None)[source] Tune the parameters of the spectral distribution. Three methods are available: log-grid-search (brute), annealing (dual_annealing) and Bayesian stochastic optimization (bso) from GpyOpt. The default object function to optimize is athena.utils.average_rrmse, which uses a cross-validation procedure from athena.utils, see Example and tutorial 06_kernel-based_AS.
- Example
>>> from athena.kas import KernelActiveSubspaces >>> from athena.feature_map import FeatureMap >>> from athena.projection_factory import ProjectionFactory >>> from athena.utils import CrossValidation, average_rrmse >>> input_dim, output_dim, n_samples = 2, 1, 30 >>> n_features, n_params = 10, 1 >>> xx = np.ones((n_samples, input_dim)) >>> f = np.ones((n_samples, output_dim)) >>> df = np.ones((n_samples, output_dim, input_dim)) >>> fm = FeatureMap(distr='laplace', bias=np.random.uniform(0, 2 * np.pi, n_features), input_dim=input_dim, n_features=n_features, params=np.zeros(n_params), sigma_f=f.var()) >>> kss = KernelActiveSubspaces(feature_map=fm, dim=1, n_features=n_features) >>> csv = CrossValidation(inputs=xx, outputs=f.reshape(-1, 1), gradients=df.reshape(n_samples, output_dim, input_dim), folds=3, subspace=kss) >>> best = fm.tune_pr_matrix(func=average_rrmse, bounds=[slice(-2, 1, 0.2) for i in range(n_params)], args=(csv, ), method='bso', maxiter=20, save_file=False)
- Parameters
func (callable) – the objective function to be minimized. Must be in the form f(x, *args), where x is the argument in the form of a 1-D array and args is a tuple of any additional fixed parameters needed to completely specify the function.
bounds (tuple) – each component of the bounds tuple must be a slice tuple of the form (low, high, step). It defines bounds for the objective function parameter in a logarithmic scale. Step will be ignored for ‘dual_annealing’ method.
args (tuple) – any additional fixed parameters needed to completely specify the objective function.
method (str) – method used to optimize the objective function. Possible values are ‘brute’, or ‘dual_annealing’. Default value is None, and the choice is made automatically wrt the dimension of self.params. If the dimension is less than 4 brute force is used, otherwise a traditional Generalized Simulated Annealing will be performed with no local search strategy applied.
maxiter (int) – the maximum number of global search iterations. Default value is 50.
save_file (bool) – True to save the optimal projection matrix
fn_args (dict) – dictionary of arguments passed to func.
- Raises
ValueError
- Returns
list that records the best score and the best projection matrix. The initial values are 0.8 and a n_features-by-input_dim numpy.ndarray of zeros.
- Return type
-
feature_map.
rff_map
(pr_matrix, bias, n_features, sigma_f)¶ Random Fourier Features map. It can be vectorized for inputs of shape n_samples-by-input_dim.
- Parameters
inputs (numpy.ndarray) – input_dim dimensional inputs to project to the RKHS.
pr_matrix (numpy.ndarray) – n_features-by-input_dim projection matrix, whose rows are sampled from the spectral distribution.
bias (numpy.ndarray) – n_features dimensional bias. It is sampled from a Unifrom distribution on the interval [0, 2*PI].
n_features (int) – dimension of the RKHS
sigma_f (int) – multiplicative term representing the empirical variance the outptus.
- Returns
n_features dimensional projection of the inputs to the RKHS
- Return type
-
feature_map.
rff_jac
(pr_matrix, bias, n_features, sigma_f)¶ Random Fourier Features map’s Jacobian. It can be vectorized for inputs of shape n_samples-by-input_dim.
- Parameters
inputs (numpy.ndarray) – input_dim dimensional inputs to project to the RKHS.
pr_matrix (numpy.ndarray) – n_features-by-input_dim projection matrix, whose rows are sampled from the spectral distribution.
bias (numpy.ndarray) – n_features dimensional bias. It is sampled from a Unifrom distribution on the interval [0, 2*PI].
n_features (int) – dimension of the RKHS
sigma_f (int) – multiplicative term representing the empirical variance the outptus.
- Returns
n_features-by-input_dim dimensional projection of the inputs to the RKHS
- Return type