skbel.tmaps

This package contains the implementation of the transport map method.

skbel.tmaps.transport

class skbel.tmaps.transport_map.TransportMap(monotone: list, nonmonotone: list, X: array, polynomial_type: str = 'hermite function', monotonicity: str = 'integrated rectifier', standardize_samples: bool = True, workers: int = 1, ST_scale_factor: float = 1.0, ST_scale_mode: str = 'dynamic', coeffs_init: float = 0.1, linearization: float = None, linearization_specified_as_quantiles: bool = True, linearization_increment: float = 1e-06, regularization: str = None, regularization_lambda: float = 0.1, quadrature_input: dict = None, rectifier_type: str = 'exponential', delta: float = 0.0)[source]

Bases: object

class Rectifier(mode: str = 'softplus', delta: float = 1e-08)[source]

Bases: object

__init__(mode: str = 'softplus', delta: float = 1e-08)[source]

This object specifies what function is used to rectify the monotone map component functions if monotonicity = ‘integrated rectifier’, before the rectifier’s output is integrated to yield a monotone map component in x_k.

Parameters:
  • mode – keyword string defining which function is used to rectify the map component functions.

  • delta – a small offset value to prevent arithmetic underflow in some rectifier functions.

Returns:

result: the rectified function value.

evaluate(X: array) array[source]

This function evaluates the specified rectifier.

Parameters:

X – an array of function evaluates to be rectified.

evaluate_dfdc(f, dfdc)[source]

This function evaluates terms used in the optimization of the map components if monotonicity = ‘separable monotonicity’.

Parameters:
  • f – an array of function evaluates to be rectified.

  • dfdc – an array of function evaluates to be rectified.

Returns:

result: the derivative of the rectified function value.

evaluate_dx(X: array) array[source]

This function evaluates the derivative of the specified rectifier.

Parameters:

X – an array of function evaluates to be rectified.

Returns:

result: the derivative of the rectified function value.

inverse(X)[source]

This function evaluates the inverse of the specified rectifier.

Parameters:

X – an array of function evaluates to be rectified.

Returns:

result: the inverse of the rectified function value.

logevaluate(X: array) array[source]

This function evaluates the logarithm of the specified rectifier.

Parameters:

X – an array of function evaluates to be rectified.

Returns:

result: the logarithm of the rectified function value.

__init__(monotone: list, nonmonotone: list, X: array, polynomial_type: str = 'hermite function', monotonicity: str = 'integrated rectifier', standardize_samples: bool = True, workers: int = 1, ST_scale_factor: float = 1.0, ST_scale_mode: str = 'dynamic', coeffs_init: float = 0.1, linearization: float = None, linearization_specified_as_quantiles: bool = True, linearization_increment: float = 1e-06, regularization: str = None, regularization_lambda: float = 0.1, quadrature_input: dict = None, rectifier_type: str = 'exponential', delta: float = 0.0)[source]

This toolbox contains functions required to construct, optimize, and evaluate transport methods.

Maximilian Ramgraber, January 2022

Parameters:
  • monotone – list specifying the structure of the monotone part of the transport map component functions.

  • nonmonotone – list specifying the structure of the non-monotone part of the transport map component functions.

  • X – N-by-D array of the training samples used to optimize the transport map, where N is the number of samples and D is the number of dimensions

  • polynomial_type – keyword which specifies what kinds of polynomials are used for the transport map component functions.

  • monotonicity – keyword which specifies through what method the transport map ensures monotonicity in the last dimensions. Must be ‘integrated rectifier’ or ‘separable monotonicity’.

  • standardize_samples – a True/False flag determining whether the transport map should standardize the training samples before optimization

  • workers – number of workers for parallel optimization. If set to 1, parallelized optimization is inactive.

  • ST_scale_factor – a float which scales the width of special terms used in the map components, such as ‘RBF 0’, ‘iRBF 0’, ‘LET 0’, or ‘RET 0’.

  • ST_scale_mode – keyword which defines whether the width of special term scale factors is determined based on neighbouring special terms (‘dynamic’) or fixed as ST_scale_factor (‘static’).

  • coeffs_init – value used to initialize the coefficients at the start of the map optimization.

  • linearization – float which specifies boundary values used to linearize the map components in the tails. Its role is specified by linearization_specified_as_quantiles.

  • linearization_specified_as_quantiles – flag which specifies whether the linearization thresholds are specified as quantiles (True) or absolute values (False). If True, boundaries are placed at linearization and 1-linearization, if False, at linearization and -linearization. Only used if linearization is not None.

  • linearization_increment – increment used for the linearization of the map component functions. Only used if linearization is not None.

  • regularization – keyword which specifies if regularization is used, and if so, what kind of regularization (‘L1’ or ‘L2’).

  • regularization_lambda – float which specifies the weight for the map coefficient regularization. Only used if regularization is not None.

  • quadrature_input – dictionary for optional keywords to overwrite the default variables in the function Gauss_quadrature. Only used if monotonicity = ‘integrated rectifier’.

  • rectifier_type – keyword which specifies which function is used to rectify the monotone map components. Only used if monotonicity = ‘integrated rectifier’.

  • delta – small increment added to the rectifier to prevent numerical underflow. Only used if monotonicity = ‘integrated rectifier’.

calculate_special_term_locations()[source]

This function calculates the location and scale parameters for special terms in the transport map definition, specifically RBF (Radial Basis Functions), iRBF (Integrated Radial Basis Functions), and LET/RET (Edge Terms).

Position and scale parameters are assigned in the order they have been defined, so make sure to define a left edge term first if you want it to be on the left side.

check_for_special_terms()[source]

This function scans through the user-provided map specifications and seeks if there are any special terms (‘RBF’, ‘iRBF’, ‘LET’, ‘RET’) among the terms of the map components.

If there are, it determines how many there are, and informs the rest of the function where these special terms should be located.

check_inputs()[source]

This function runs some preliminary checks on the input provided, alerting the user to any possible input errors.

function_constructor_alternative()[source]

This function assembles the string for the monotone and non-monotone map components, then converts these strings into functions.

function_derivative_constructor_alternative()[source]

This function is the complement to ‘function_constructor_alternative’, but instead constructs the derivative of the map’s component functions.

It constructs the functions’ strings, then converts them into functions.

gauss_quadrature(f, a, b, order=100, args=None, Ws=None, xis=None, adaptive=False, threshold=1e-06, increment=1, verbose=False, full_output=False)[source]

This function implements a Gaussian quadrature numerical integration scheme. It is used if the monotonicity = ‘integrated rectifier’, for which monotonicity is ensured by integrating a strictly positive function obtained from a rectifier.

Parameters:
  • f – function to be integrated element-wise.

  • a – lower bound for integration, defined as either a scalar or a vector.

  • b – upper bound for integration, defined as either a scalar or a vector.

  • order – order of the Legendre polynomial used for the integration scheme.

  • args – a dictionary with supporting keyword arguments to be passed to the function.

  • Ws – weights of the integration points, can be calculated in advance to speed up the computation. Is calculated by the integration scheme, if not specified.

  • xis – positions of the integration points, can be calculated in advance to speed up the computation. Is calculated by the integration scheme, if not specified.

  • full_output – Flag for whether the positions and weights of the integration points should be returned along with the integration results. If True, returns a tuple with (results,order,xis,Ws). If False, only returns results.

  • adaptive – flag which determines whether the numerical scheme should adjust the order of the Legendre polynomial adaptively (True) or use the integer provided by ‘order’ (False).

  • threshold – threshold for the difference in the adaptive integration, adaptation stops after difference in integration

  • increment – increment by which the order is increased in each adaptation cycle. Higher values correspond to larger steps.

  • verbose – flag which determines whether information about the integration process should be printed to console (True) or not (False).

Returns:

the result of the integration, either as a scalar or a vector.

inverse_map(Y: array, X_precalc: array = None)[source]

This function evaluates the inverse transport map, mapping samples from a multivariate standard Gaussian back to the target distribution. If X_precalc is specified, the map instead evaluates a conditional of the target distribution given X_precalc. The function assumes any precalculated output are the FIRST dimensions of the total output. If X_precalc is specified, its dimensions and the input dimensions must sum to the full dimensionality of sample space.

Parameters:
  • Y – N-by-D or N-by-(D-E) array of reference distribution samples to be mapped to the target distribution, where N is the number of samples, D is the number of target distribution dimensions, and E the number of pre-specified dimensions (if X_precalc is specified).

  • X_precalc – N-by-E array of samples in the space of the target distribution, used to condition the lower D-E dimensions during the inversion process.

Returns:

N-by-D array of samples in the space of the target distribution.

map(X: array = None)[source]

This function maps the samples X from the target distribution to the standard multivariate Gaussian reference distribution. If X has not been provided, the samples in storage will be used instead.

Parameters:

X – N-by-D array of the training samples used to optimize the transport map, where N is the number of samples and D is the number of dimensions.

Returns:

N-by-D array of the mapped samples

objective_function(coeffs: array, k: int, div: int = 0)[source]

This function evaluates the objective function used in the optimization of the map’s component functions.

Parameters:
  • coeffs – a vector containing the coefficients for both the non-monotone and monotone terms of the k-th map component function. Is replaced for storage is specified as None.

  • k – an integer variable defining what map component

  • div – an integer specifying where the cutoff between the non-monotone and monotone coefficients in ‘coeffs’ is.

Returns:

the value of the objective function

objective_function_jacobian(coeffs: array, k: int, div: int = 0)[source]

This function evaluates the derivative of the objective function used in the optimization of the map’s component functions.

Parameters:
  • coeffs – a vector containing the coefficients for both the non-monotone and monotone terms of the k-th map component function. Is replaced for storage is specified as None.

  • k – an integer variable defining what map component is being evaluated. Corresponds to a dimension of sample space.

  • div – an integer specifying where the cutoff between the non-monotone and monotone coefficients in ‘coeffs’ is.

optimize()[source]

This function optimizes the map’s component functions, seeking the coefficients which best map the samples to a standard multivariate Gaussian distribution.

precalculate()[source]

This function pre-calculates matrices of basis function evaluations for the samples provided.

These matrices can be used to optimize the maps more quickly.

reset(X: array)[source]

This function is used if the transport map has been initiated with a different set of samples. It resets the standardization variables and the map’s coefficients, requiring new optimization.

Parameters:

X – N-by-D array of the training samples used to optimize the transport map, where N is the number of samples and D is the number of dimensions

s(x: array, k: int, coeffs_nonmon: array = None, coeffs_mon: array = None)[source]

This function evaluates the k-th map component.

Parameters:
  • x – N-by-D array of the training samples used to optimize the transport map, where N is the number of samples and D is the number of dimensions. Can be None, at which point it is replaced with X from storage.

  • k – an integer variable defining what map component is being evaluated. Corresponds to a dimension of sample space.

  • coeffs_nonmon – a vector specifying the coefficients of the non-monotone part of the map component’s terms, i.e., those entries which do not depend on x_k. This vector is replaced from storage if it is not overwritten.

  • coeffs_mon – a vector specifying the coefficients of the monotone part of the map component’s terms, i.e., those entries which do not depend on x_k. This vector is replaced from storage if it is not overwritten.

Returns:

N-by-1 array of the k-th map component evaluated at the samples in x

vectorized_root_search_bisection(X: array, Yk: array, k: int)[source]

This function searches for the roots of the k-th map component through bisection. It is called in the inverse_map function.

Parameters:
  • X – N-by-k array of samples inverted so far, where the k-th column still contains the reference samples used as a residual in the root finding process

  • Yk – a vector containing the target values in the k-th dimension, for which the root finding algorithm must solve.

  • k – an integer variable defining what map component is being evaluated. Corresponds to a dimension of sample space.

worker_task(k: int, task_supervisor: list)[source]

This function provides the optimization task for the k-th map component function to a worker (if parallelization is used), or applies it in sequence (if no parallelization is used). This specific function only becomes active if monotonicity = ‘integrated rectifier’.

Parameters:
  • k – an integer variable defining what map component is being evaluated. Corresponds to a dimension of sample space.

  • task_supervisor

    a shared list which informs the main process how many optimization tasks have already been computed. This list should not be specified by the user, it only serves to provide information about the

    optimization progress.

worker_task_monotone(k: int, task_supervisor: list)[source]

This function provides the optimization task for the k-th map component function to a worker (if parallelization is used), or applies it in sequence (if no parallelization is used). This specific function only becomes active if monotonicity = ‘separable monotonicity’.

Parameters:
  • k – an integer variable defining what map component is being evaluated. Corresponds to a dimension of sample space.

  • task_supervisor – a shared list which informs the main process how many optimization tasks have already been computed. This list should not be specified by the user, it only serves to provide information about the optimization progress.

Returns:

a list containing the optimized coefficients for the k-th map component function.

write_basis_function(term: list, mode: str = 'standard', k: int = None)[source]

This function assembles a string for a specific term of the map component functions. This can be a polynomial, a Hermite function, a radial basis function, or similar.

Parameters:
  • term – either an empty list, a list, or a string, which specifies what function type this term is supposed to be. Empty lists correspond to a constant, lists specify polynomials or Hermite functions, and strings denote special terms such as radial basis functions.

  • mode – a keyword which defines whether the term’s string returned should be conventional (‘standard’) or its derivative (‘derivative’).

  • k – an integer specifying what dimension of the samples the ‘term’ corresponds to. Used to clarify with respect to what variable we take the derivative.