Module API

Module arby.rom

Reduced Order Modeling module.

class arby.rom.ReducedOrderModel(training_set, physical_points, parameter_points, integration_rule: str = 'riemann', greedy_tol: float = 1e-12, poly_deg: int = 3)

Build reduced order models from training data.

This class comprises a set of tools to build and handle reduced bases, empirical interpolants and predictive models from pre-computed training set of functions. The underlying or ground truth model describing the training set is a real function g(v,x) parameterized by a training parameter v. The physical variable x belongs to a domain for which an inner product can defined. The surrogate model is built bringing together the Reduced Basis (RB) greedy algorithm and the Empirical Interpolation Method (EIM) to work in synergy towards a predictive model for the ground truth model.

Parameters
  • training_set (array_like) – Array of training functions.

  • physical_points (array_like) – Array of physical points.

  • parameter_points (array_like) – Array of parameter points.

  • integration_rule (str, optional) – The quadrature rule to define an integration scheme. Default = “riemann”.

  • greedy_tol (float, optional) – The greedy tolerance as a stopping condition for the reduced basis greedy algorithm. Default = 1e-12.

  • poly_deg (int, optional) – Degree <= 5 of polynomials used for splines. Default = 3.

Examples

Build a surrogate model

>>> from arby import ReducedOrderModel as ROM

Input the three most important parameters (the others are optional).

>>> model = ROM(training_set, physical_points, parameter_points)

Build/evaluate the surrogate model. The building stage is done once and for all at the first call. It could take some time for large training sets. For this reason it is called the offline stage. Subsequent calls will invoke the already built surrogate model and just evaluates it. That corresponds to the online stage.

>>> model.surrogate(parameter)

For attempting to improve the model’s accuracy without the addition of more training functions, tune the class parameters greedy_tol and poly_deg to control the precision of the reduced basis (see the arby.reduced_basis method) or the internal splines model.

property Ntrain_

Return the number of training functions or parameter points.

property Nsamples_

Return the number of samples or physical points.

property basis_

Return a reduced basis object.

The reduced basis is computed at the first call and stored as a class object of arby.Basis, which comprises several tools for handling bases. See also the arby.reduced_basis documentation.

property greedy_indices_

Greedy indices from the RB algorithm.

See the arby.reduced_basis documentation.

property greedy_errors_

Errors computed in the RB algorithm.

See the arby.reduced_basis documentation.

property projection_matrix_

Projection coefficients from the RB algorithm.

See the arby.reduced_basis documentation.

property eim_

Return EIM data.

See arby.Basis.eim_ documentation.

surrogate(param)

Evaluate the surrogate model at parameter/s.

Build a surrogate model valid for the entire parameter domain. The building stage is performed only once for the first function call. For subsequent calls, the method invokes the already fitted model and just evaluates it. The output is an array storing surrogate evaluations at the parameter/s.

Parameters

param (float or array_like(float)) – Point or set of parameters inside the parameter domain.

Returns

h_surrogate – The evaluated surrogate function for the given parameters.

Return type

numpy.ndarray

Module arby.basis

Basis analysis module.

class arby.basis.EIM(interpolant: numpy.ndarray, nodes: list)

Container for EIM information.

Parameters
  • interpolant (numpy.ndarray) – Interpolant matrix.

  • nodes (list) – EIM nodes.

class arby.basis.RB(basis: numpy.ndarray, indices: numpy.ndarray, errors: numpy.ndarray, projection_matrix: numpy.ndarray)

Container for RB information.

Parameters
  • basis (np.ndarray) – Reduced basis object.

  • indices (np.ndarray) – Greedy indices.

  • errors (np.ndarray) – Greedy errors.

  • projection_matrix (np.ndarray) – Projection coefficients.

class arby.basis.Basis(data, integration: numpy.ndarray)

Basis object and utilities.

Create a basis object introducing an orthonormalized set of functions data and an integration class instance to enable integration utilities for the basis.

Parameters
property Nbasis_: int

Return the number of basis elements.

property eim_: arby.basis.EIM

Implement EIM algorithm.

The Empirical Interpolation Method (EIM) [TiglioAndVillanueva2021] introspects the basis and selects a set of interpolation nodes from the physical domain for building an interpolant matrix using the basis and the selected nodes. The interpolant matrix can be used to approximate a field of functions for which the span of the basis is a good approximant.

Returns

Container for EIM data. Contains (interpolant, nodes).

Return type

arby.basis.EIM

projection_error(h, s=(None))

Compute the squared projection error of a function h onto the basis.

The error is computed in the L2 norm (continuous case) or the 2-norm (discrete case), that is, ||h - P h||^2, where P denotes the projector operator associated to the basis.

Parameters
  • h (np.ndarray) – Function to be projected.

  • s (tuple, optional) – Slice the basis. If the slice is not provided, the whole basis is considered. Default = (None,)

Returns

error – Square of the projection error.

Return type

float

project(h, s=(None))

Project a function h onto the basis.

This method represents the action of projecting the function h onto the span of the basis.

Parameters
  • h (np.ndarray) – Function or set of functions to be projected.

  • s (tuple, optional) – Slice the basis. If the slice is not provided, the whole basis is considered. Default = (None,)

Returns

projected_function – Projection of h onto the basis.

Return type

np.ndarray

interpolate(h)

Interpolate a function h at EIM nodes.

This method uses the basis and associated EIM nodes (see the arby.Basis.eim_ method) for interpolation.

Parameters

h (np.ndarray) – Function or set of functions to be interpolated.

Returns

h_interpolated – Interpolated function at EIM nodes.

Return type

np.ndarray

arby.basis.gram_schmidt(functions, integration, max_iter=3)numpy.ndarray

Orthonormalize a set of functions.

This algorithm implements the Iterated, Modified Gram-Schmidt (GS) algorithm [Hoffmann1989] to build an orthonormal basis from a set of functions.

Parameters
  • functions (array_like, shape=(m, L)) – Functions to be orthonormalized, where m is the number of functions and L is the length of the sample.

  • integration (arby.integrals.Integration) – Instance of the Integration class for defining inner products.

  • max_iter (int, optional) – Maximum number of iterations. Default = 3.

Returns

basis – Orthonormalized array.

Return type

numpy.ndarray

Raises

ValueError – If functions are not linearly independent for a given tolerance.

References

Hoffmann1989

Hoffmann, W. Iterative algorithms for Gram-Schmidt orthogonalization. Computing 41, 335-348 (1989). https://doi.org/10.1007/BF02241222

arby.basis.reduced_basis(training_set, physical_points, integration_rule='riemann', greedy_tol=1e-12, normalize=False)arby.basis.RB

Build a reduced basis from training data.

This function implements the Reduced Basis (RB) greedy algorithm for building an orthonormalized reduced basis out from training data. The basis is built for reproducing the training functions within a user specified tolerance [TiglioAndVillanueva2021] by linear combinations of its elements. Tuning the greedy_tol parameter allows to control the representation accuracy of the basis.

The integration_rule parameter specifies the rule for defining inner products. If the training functions (rows of the training_set) does not correspond to continuous data (e.g. time), choose "euclidean". Otherwise choose any of the quadratures defined in the arby.Integration class.

Set the boolean normalize to True if you want to normalize the training set before running the greedy algorithm. This condition not only emphasizes on structure over scale but may leads to noticeable speedups for large datasets.

The output is a container which comprises RB data: a basis object storing the reduced basis and handling tools (see arby.Basis); the greedy errors corresponding to the maxima over the training set of the squared projection errors for each greedy swept; the greedy indices locating the most relevant training functions used for building the basis; and the projection_matrix storing the projection coefficients generated by the greedy algorithm. For example, we can recover the training set (more precisely, a compressed version of it) by multiplying the projection matrix with the reduced basis.

Parameters
  • training_set (numpy.ndarray) – The training set of functions.

  • physical_points (numpy.ndarray) – Physical points for quadrature rules.

  • integration_rule (str, optional) – The quadrature rule to define an integration scheme. Default = “riemann”.

  • greedy_tol (float, optional) – The greedy tolerance as a stopping condition for the reduced basis greedy algorithm. Default = 1e-12.

  • normalize (bool, optional) – True if you want to normalize the training set. Default = False.

Returns

Container for RB data. Contains (basis, errors, indices, projection_matrix).

Return type

arby.basis.RB

Raises

ValueError – If training_set.shape[1] doesn’t coincide with quadrature rule weights.

Notes

If normalize is True, the projection coefficients are with respect to the original basis but the greedy errors are relative to the normalized training set.

References

TiglioAndVillanueva2021(1,2)

Reduced Order and Surrogate Models for Gravitational Waves. Tiglio, M. and Villanueva A. arXiv:2101.11608 (2021)

Module arby.integrals

Integration schemes module.

class arby.integrals.Integration(interval, rule='riemann')

Integration scheme.

This class fixes a frame for performing integrals, inner products and derived operations. An integral is defined by a quadrature rule composed by nodes and weights which are used to construct a discrete approximation to the true integral (or inner product).

For completeness, an “euclidean” rule is available for which inner products reduce to simple discrete dot products.

Parameters
  • interval (numpy.ndarray) – Equispaced set of points as domain for integrals or inner products.

  • rule (str, optional) – Quadrature rule. Default = “riemann”. Available = (“riemann”, “trapezoidal”, “euclidean”)

integral(f)

Integrate a function.

Parameters

f (np.ndarray) – Real or complex numbers array.

dot(f, g)

Return the dot product between functions.

Parameters
  • f (np.ndarray) – Real or complex numbers array.

  • g (np.ndarray) – Real or complex numbers array.

norm(f)

Return the norm of a function.

Parameters

f (np.ndarray) – Real or complex numbers array.

normalize(f)

Normalize a function.

Parameters

f (np.ndarray) – Real or complex numbers array.