TransportMaps.Functionals¶
Classes

Class | Description |
---|---|
Function | Abstract class for parametric approximation \(f_{\bf a}:\mathbb{R}^d\rightarrow\mathbb{R}\) of \(f:\mathbb{R}^d\rightarrow\mathbb{R}\). |
ParametricFunctionApproximation | Abstract class for parametric approximation \(f_{\bf a}:\mathbb{R}^d\rightarrow\mathbb{R}\) of \(f:\mathbb{R}^d\rightarrow\mathbb{R}\). |
TensorizedFunctionApproximation | [Abstract] Class for approximations using tensorization of unidimensional basis |
LinearSpanApproximation | Parametric function \(f_{\bf a} = \sum_{{\bf i} \in \mathcal{I}} {\bf a}_{\bf i} \Phi_{\bf i}\) |
IntegratedSquaredParametricFunctionApproximation | Parameteric function \(f_{\bf a}({\bf x}) = \int_0^{x_d} h_{\bf a}^2(x_1,\ldots,x_{d-1},t) dt\) |
MonotonicFunctionApproximation | Abstract class for the approximation \(f \approx f_{\bf a} = \sum_{{\bf i} \in \mathcal{I}} {\bf a}_{\bf i} \Phi_{\bf i}\) assumed to be monotonic in \(x_d\) |
MonotonicLinearSpanApproximation | Approximation of the type \(f \approx f_{\bf a} = \sum_{{\bf i} \in \mathcal{I}} {\bf a}_{\bf i} \Phi_{\bf i}\), monotonic in \(x_d\) |
MonotonicIntegratedExponentialApproximation | Integrated Exponential approximation. |
MonotonicIntegratedSquaredApproximation | Integrated Squared approximation. |
ProductDistributionParametricPullbackComponentFunction | Parametric function \(f[{\bf a}](x_{1:k}) = \log\pi\circ T_k[{\bf a}](x_{1:k}) + \log\partial_{x_k}T_k[{\bf a}](x_{1:k})\) |
MonotonicFrozenFunction | [Abstract] Frozen function. No optimization over the coefficients allowed. |
FrozenLinear | Frozen Linear map \({\bf x} \rightarrow a_1 + a_2 {\bf x}_d\) |
FrozenExponential | Frozen Exponential map \(f_{\bf a}:{\bf x} \mapsto \exp( {\bf x}_d )\) |
FrozenGaussianToUniform | Frozen Gaussian To Uniform map. |
Documentation
-
class
TransportMaps.Functionals.
Function
(dim)[source]¶ Abstract class for parametric approximation \(f_{\bf a}:\mathbb{R}^d\rightarrow\mathbb{R}\) of \(f:\mathbb{R}^d\rightarrow\mathbb{R}\).
Parameters: dim (int) – number of dimensions -
evaluate
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: (
ndarray
[\(m\)]) – function evaluations- x (
-
grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
- x (
-
grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
partial2_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\partial^2_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial^2_{x_d} f_{\bf a}({\bf x})\)
- x (
-
partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
-
class
TransportMaps.Functionals.
ParametricFunctionApproximation
(dim)[source]¶ Abstract class for parametric approximation \(f_{\bf a}:\mathbb{R}^d\rightarrow\mathbb{R}\) of \(f:\mathbb{R}^d\rightarrow\mathbb{R}\).
Parameters: dim (int) – number of dimensions -
evaluate
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: (
ndarray
[\(m\)]) – function evaluations- x (
-
grad_a
(x, precomp=None, idxs_slice=slice(None, None, None))[source]¶ [Abstract] Evaluate \(\nabla_{\bf a} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a} f_{\bf a}({\bf x})\)
- x (
-
grad_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
hess_a
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla^2_{\bf a} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a} f_{\bf a}({\bf x})\)
- x (
-
hess_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
n_coeffs
¶ [Abstract] Get the number \(N\) of coefficients \({\bf a}\)
Returns: ( int
) – number of coefficients
-
partial2_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\partial^2_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial^2_{x_d} f_{\bf a}({\bf x})\)
- x (
-
partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ [Abstract] Evaluate \(\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
precomp_evaluate
(x, precomp=None)[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_grad_x
(x, precomp=None)[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_grad_x_partial_xd
(x, precomp=None)[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_partial2_xd
(x, precomp=None)[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\partial^2_{x_d} f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_partial_xd
(x, precomp=None)[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\partial_{x_d} f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
regression_action_storage_hess_a_objective
(a, v, params)[source]¶ Assemble/fetch Hessian \(\nabla_{\bf a}^2 \Vert f - f_{\bf a} \Vert^2_{\pi}\) and evaluate its action on \(v\)
Parameters:
-
regression_grad_a_objective
(a, params)[source]¶ Objective function \(\nabla_{\bf a} \Vert f - f_{\bf a} \Vert^2_{\pi}\)
Parameters:
-
-
class
TransportMaps.Functionals.
TensorizedFunctionApproximation
(basis_list, full_basis_list=None)[source]¶ [Abstract] Class for approximations using tensorization of unidimensional basis
Parameters: -
precomp_evaluate
(x, precomp=None, precomp_type='uni')[source]¶ Precompute the uni-variate Vandermonde matrices for the evaluation of \(f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_grad_x
(x, precomp=None)[source]¶ Precompute the uni-variate Vandermonde matrices for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Letting \(\Phi^{(i)}(x_i)\) being the uni-variate Vandermonde in \(x_i\), the
i
-th element of the returned list is \(\partial_{x_i}\Phi^{(i)}(x_i)\).Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_grad_x_partial_xd
(x, precomp=None)[source]¶ Precompute uni-variate Vandermonde matrices for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_hess_x
(x, precomp=None)[source]¶ Precompute the uni-variate Vandermonde matrices for the evaluation of \(\nabla^2_{\bf x} f_{\bf a}\) at
x
Letting \(\Phi^{(i)}(x_i)\) being the uni-variate Vandermonde in \(x_i\), the
i
-th element of the returned list is \(\partial^2_{x_i}\Phi^{(i)}(x_i)\).Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_hess_x_partial_xd
(x, precomp=None)[source]¶ Precompute uni-variate Vandermonde matrices for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_partial2_xd
(x, precomp=None)[source]¶ Precompute uni-variate Vandermonde matrix for the evaluation of \(\partial^2_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_partial3_xd
(x, precomp=None)[source]¶ Precompute uni-variate Vandermonde matrix for the evaluation of \(\partial^3_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
-
class
TransportMaps.Functionals.
LinearSpanApproximation
(basis_list, spantype=None, order_list=None, multi_idxs=None, full_basis_list=None)[source]¶ Parametric function \(f_{\bf a} = \sum_{{\bf i} \in \mathcal{I}} {\bf a}_{\bf i} \Phi_{\bf i}\)
Parameters: - basis_list (list) – list of \(d\)
OrthogonalBasis
- spantype (str) – Span type. ‘total’ total order, ‘full’ full order, ‘midx’ multi-indeces specified
- order_list (
list
ofint
) – list of orders \(\{N_i\}_{i=0}^d\) - multi_idxs (list) – list of tuples containing the active multi-indices
- full_basis_list (list) – full list of
Basis
.basis_list
is a subset offull_basis_list
. This may be used to automatically increase the input dimension of the approximation.
-
evaluate
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: (
ndarray
[\(m\)]) – function evaluations- x (
-
get_directional_orders
()[source]¶ Get the maximum orders of the univariate part of the representation.
Returns: ( list
[d]int
) – orders
-
grad_a
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a} f_{\bf a}({\bf x})\)
- x (
-
grad_a_grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a} \nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N,d\)]) – \(\nabla_{\bf a} \nabla_{\bf x} f_{\bf a}({\bf x})\)
- x (
-
grad_a_grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a}\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N,d\)]) – \(\nabla_{\bf a}\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
grad_a_hess_x
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a} \nabla^2_{\bf x} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N,d,d\)]) – \(\nabla_{\bf a} \nabla^2_{\bf x} f_{\bf a}({\bf x})\)
- x (
-
grad_a_hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwars)[source]¶ Evaluate \(\nabla_{\bf a}\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla_{\bf a}\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
grad_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
grad_a_squared
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\left(\nabla_{\bf a} f_{\bf a}\right)\left(\nabla_{\bf a} f_{\bf a}\right)^T\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N,N\)]) – \(\left(\nabla_{\bf a} f_{\bf a}\right)\left(\nabla_{\bf a} f_{\bf a}\right)^T\)
- x (
-
grad_a_t
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\left(\nabla_{\bf a} f_{\bf a}\right)^T\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\left(\nabla_{\bf a} f_{\bf a}({\bf x}\right)^T)\)
- x (
-
grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
- x (
-
grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
hess_a
(x, precomp=None, idxs_slice=slice(None, None, None), *arg, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf a} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a} f_{\bf a}({\bf x})\)
- (
-
hess_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
hess_x
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla^2_{\bf x} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x} f_{\bf a}({\bf x})\)
- x (
-
hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
partial2_xd
(x, precomp=None, idxs_slice=slice(None, None, None))[source]¶ Evaluate \(\partial^2_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m\)]) – \(\partial^2_{x_d} f_{\bf a}({\bf x})\)
- (
-
partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
precomp_VVT_evaluate
(x, precomp=None)[source]¶ Precompute the product \(VV^T\) of the multi-variate Vandermonde matrices for the evaluation of \(f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_Vandermonde_evaluate
(x, precomp=None)[source]¶ Precompute the multi-variate Vandermonde matrices for the evaluation of \(f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_Vandermonde_grad_x
(x, precomp=None)[source]¶ Precompute the multi-variate Vandermonde matrices for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns:
-
precomp_Vandermonde_grad_x_partial_xd
(x, precomp=None)[source]¶ Precompute multi-variate Vandermonde matrices for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_Vandermonde_hess_x
(x, precomp=None)[source]¶ Precompute the multi-variate Vandermonde matrices for the evaluation of \(\nabla^2_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_Vandermonde_hess_x_partial_xd
(x, precomp=None)[source]¶ Precompute Vandermonde matrices for the evaluation of \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_Vandermonde_partial2_xd
(x, precomp=None)[source]¶ Precompute multi-variate Vandermonde matrix for the evaluation of \(\partial^2_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_Vandermonde_partial_xd
(x, precomp=None)[source]¶ Precompute multi-variate Vandermonde matrix for the evaluation of \(\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns:
-
precomp_regression
(x, precomp=None, *args, **kwargs)[source]¶ Precompute necessary structures for the speed up of
regression()
Parameters: Returns: (
dict
) – dictionary of necessary strucutres
-
regression
(f, fparams=None, d=None, qtype=None, qparams=None, x=None, w=None, x0=None, regularization=None, tol=0.0001, maxit=100, batch_size=(None, None), mpi_pool=None, import_set=set())[source]¶ Compute \({\bf a}^* = \arg\min_{\bf a} \Vert f - f_{\bf a} \Vert_{\pi}\).
Parameters: - f (
Function
orndarray
[\(m\)]) – function \(f\) or its functions values - d (Distribution) – distribution \(\pi\)
- fparams (dict) – parameters for function \(f\)
- qtype (int) – quadrature type to be used for the approximation of \(\mathbb{E}_{\pi}\)
- qparams (object) – parameters necessary for the construction of the quadrature
- x (
ndarray
[\(m,d\)]) – quadrature points used for the approximation of \(\mathbb{E}_{\pi}\) - w (
ndarray
[\(m\)]) – quadrature weights used for the approximation of \(\mathbb{E}_{\pi}\) - x0 (
ndarray
[\(N\)]) – coefficients to be used as initial values for the optimization - regularization (dict) – defines the regularization to be used.
If
None
, no regularization is applied. If keytype=='L2'
then applies Tikonhov regularization with coefficient in keyalpha
. - tol (float) – tolerance to be used to solve the regression problem.
- maxit (int) – maximum number of iterations
- batch_size (
list
[2] ofint
) – the list contains the size of the batch to be used for each iteration. A size1
correspond to a completely non-vectorized evaluation. A sizeNone
correspond to a completely vectorized one. - mpi_pool (
mpi_map.MPI_Pool
) – pool of processes to be used - import_set (set) – list of couples
(module_name,as_field)
to be imported asimport module_name as as_field
(for MPI purposes)
Returns: (
tuple
[\(N\)],list
)) – containing the \(N\) coefficients and log information from the optimizer.See also
TransportMaps.TriangularTransportMap.regression()
Note
the resulting coefficients \({\bf a}\) are automatically set at the end of the optimization. Use
coeffs()
in order to retrieve them.Note
The parameters
(qtype,qparams)
and(x,w)
are mutually exclusive, but one pair of them is necessary.- f (
- basis_list (list) – list of \(d\)
-
class
TransportMaps.Functionals.
IntegratedSquaredParametricFunctionApproximation
(h, integ_ord_mult=6)[source]¶ Parameteric function \(f_{\bf a}({\bf x}) = \int_0^{x_d} h_{\bf a}^2(x_1,\ldots,x_{d-1},t) dt\)
Parameters: - h (
ParametricFunctionApproximation
) – parametric function \(h\) - integ_ord_mult (int) – multiplier for the number of Gauss points to be used
in the approximation of \(\int_0^{{\bf x}_d}\). The resulting number of
points is given by the product of the order in the \(d\) direction
and
integ_ord_mult
.
-
evaluate
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: (
ndarray
[\(m\)]) – function evaluations- x (
-
grad_a
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a} f_{\bf a}({\bf x})\)
- x (
-
grad_a_grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla{\bf a} \nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d\)]) – \(\nabla_{\bf a}\nabla_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_a_grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a}\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m,N,d\)]) – \(\nabla_{\bf a}\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
grad_a_hess_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla{\bf a} \nabla^2_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d,d\)]) – \(\nabla{\bf a} \nabla^2_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_a_hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a}\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m,N,d,d\)]) – \(\nabla_{\bf a}\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
grad_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
hess_a
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf a} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a} f_{\bf a}({\bf x})\)
- (
-
hess_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
hess_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x} f_{\bf a}({\bf x})\)
- (
-
hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
partial2_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\partial^2_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial^2_{x_d} f_{\bf a}({\bf x})\)
- x (
-
partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
precomp_evaluate
(x, precomp=None, precomp_type='uni')[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_grad_x
(x, precomp, precomp_type='uni')[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_grad_x_partial_xd
(x, precomp=None, precomp_type='uni')[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_partial2_xd
(x, precomp=None, precomp_type='uni')[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\partial^2_{x_d} f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
- h (
-
class
TransportMaps.Functionals.
MonotonicFunctionApproximation
(dim)[source]¶ Abstract class for the approximation \(f \approx f_{\bf a} = \sum_{{\bf i} \in \mathcal{I}} {\bf a}_{\bf i} \Phi_{\bf i}\) assumed to be monotonic in \(x_d\)
The class defines a series of methods peculiar to monotonic functions.
-
allocate_cache_minimize_kl_divergence_component
(x)[source]¶ Allocate cache space for the KL-divergence minimization
Parameters: x ( ndarray
[\(m,d\)]) – evaluation points
-
inverse
(xmd, y, xtol=1e-12, rtol=1e-15)[source]¶ Compute \({\bf x}_d\) s.t. \(f_{\bf a}({\bf x}_{1:d-1},{\bf x}_d) - y = 0\).
Given the fixed coordinates \({\bf x}_{1:d-1}\), the value \(y\), find the last coordinate \({\bf x}_d\) such that:
\[f_{\bf a}({\bf x}_{1:d-1},{\bf x}_d) - y = 0\]We will define this value the inverse of \(f_{\bf a}({\bf x})\) and denote it by \(f_{\bf a}^{-1}({\bf x}_{1:d-1})(y)\).
Parameters: Returns: (
float
) – inverse value \(x\).
-
minimize_kl_divergence_component
(f, x, w, x0=None, regularization=None, tol=0.0001, maxit=100, ders=2, fungrad=False, precomp_type='uni', batch_size=None, cache_level=1, mpi_pool=None)[source]¶ Compute \({\bf a}^\star = \arg\min_{\bf a}-\sum_{i=0}^m \log\pi\circ T_k(x_i) + \log\partial_{x_k}T_k(x_i) = \arg\min_{\bf a}-\sum_{i=0}^m f(x_i)\)
Parameters: - f (ProductDistributionParametricPullbackComponentFunction) – function \(f\)
- x (
ndarray
[\(m,d\)]) – quadrature points - w (
ndarray
[\(m\)]) – quadrature weights - x0 (
ndarray
[\(N\)]) – coefficients to be used as initial values for the optimization - regularization (dict) – defines the regularization to be used.
If
None
, no regularization is applied. If keytype=='L2'
then applies Tikonhov regularization with coefficient in keyalpha
. - tol (float) – tolerance to be used to solve the KL-divergence problem.
- maxit (int) – maximum number of iterations
- ders (int) – order of derivatives available for the solution of the optimization problem. 0 -> derivative free, 1 -> gradient, 2 -> hessian.
- fungrad (bool) – whether the distributions \(\pi_1,\pi_2\) provide the method
Distribution.tuple_grad_x_log_pdf()
computing the evaluation and the gradient in one step. This is used only forders==1
. - precomp_type (str) – whether to precompute univariate Vandermonde matrices ‘uni’ or multivariate Vandermonde matrices ‘multi’
- batch_size (
list
[3 or 2] ofint
orlist
ofbatch_size
) – the list contains the size of the batch to be used for each iteration. A size1
correspond to a completely non-vectorized evaluation. A sizeNone
correspond to a completely vectorized one. If the target distribution is aProductDistribution
, then the optimization problem decouples andbatch_size
is a list of lists containing the batch sizes to be used for each component of the map. - cache_level (int) – use high-level caching during the optimization, storing the
function evaluation
0
, and the gradient evaluation1
or nothing-1
- mpi_pool (
mpi_map.MPI_Pool
orlist
ofmpi_pool
) – pool of processes to be used,None
stands for one process. If the target distribution is aProductDistribution
, then the minimization problem decouples andmpi_pool
is a list containing ``mpi_pool``s for each component of the map.
-
minimize_kl_divergence_component_grad_a_objective
(a, params)[source]¶ Gradient of the objective function \(-\sum_{i=0}^m \nabla_{\bf a} f[{\bf a}](x_i) = -\sum_{i=0}^m \nabla_{\bf a} \left( \log\pi\circ T_k[{\bf a}](x_i) + \log\partial_{x_k}T_k[{\bf a}](x_i)\right)\)
Parameters:
-
minimize_kl_divergence_component_hess_a_objective
(a, params)[source]¶ Hessian of the objective function \(-\sum_{i=0}^m \nabla^2_{\bf a} f[{\bf a}](x_i) = -\sum_{i=0}^m \nabla^2_{\bf a} \left( \log\pi\circ T_k[{\bf a}](x_i) + \log\partial_{x_k}T_k[{\bf a}](x_i)\right)\)
Parameters:
-
minimize_kl_divergence_component_objective
(a, params)[source]¶ Objective function \(-\sum_{i=0}^m f(x_i) = -\sum_{i=0}^m \log\pi\circ T_k(x_i) + \log\partial_{x_k}T_k(x_i)\)
Parameters:
-
partial_xd_inverse
(xmd, y)[source]¶ Compute \(\partial_y f_{\bf a}^{-1}({\bf x}_{1:d-1})(y)\).
Parameters: Returns: (
float
) – derivative of the inverse value \(x\).
-
partial_xd_misfit
(x, args)[source]¶ Compute \(\partial_{x_d} f_{\bf a}({\bf x}) - y = \partial_{x_d} f_{\bf a}({\bf x})\)
Given the fixed coordinates \({\bf x}_{1:d-1}\), the value \(y\), and the last coordinate \({\bf x}_d\), compute:
\[\partial f_{\bf a}({\bf x}_{1:d-1},{\bf x}_d)\]Parameters: Returns: (
float
) – misfit derivative.
-
precomp_minimize_kl_divergence_component
(x, params, precomp_type='uni')[source]¶ Precompute necessary structures for the speed up of
minimize_kl_divergence_component()
Parameters: Returns: - (
tuple
(None,:class:dict<dict>)) – dictionary of necessary strucutres. The first argument is needed for consistency with
- (
-
regression
(f, fparams=None, d=None, qtype=None, qparams=None, x=None, w=None, x0=None, regularization=None, tol=0.0001, maxit=100, batch_size=(None, None, None), mpi_pool=None, import_set=set())[source]¶ Compute \({\bf a}^* = \arg\min_{\bf a} \Vert f - f_{\bf a} \Vert_{\pi}\).
Parameters: - f (
Function
orndarray
[\(m\)]) – function \(f\) or its functions values - fparams (dict) – parameters for function \(f\)
- d (Distribution) – distribution \(\pi\)
- qtype (int) – quadrature type to be used for the approximation of \(\mathbb{E}_{\pi}\)
- qparams (object) – parameters necessary for the construction of the quadrature
- x (
ndarray
[\(m,d\)]) – quadrature points used for the approximation of \(\mathbb{E}_{\pi}\) - w (
ndarray
[\(m\)]) – quadrature weights used for the approximation of \(\mathbb{E}_{\pi}\) - x0 (
ndarray
[\(N\)]) – coefficients to be used as initial values for the optimization - regularization (dict) – defines the regularization to be used.
If
None
, no regularization is applied. If keytype=='L2'
then applies Tikonhov regularization with coefficient in keyalpha
. - tol (float) – tolerance to be used to solve the regression problem.
- maxit (int) – maximum number of iterations
- batch_size (
list
[3] ofint
) – the list contains the size of the batch to be used for each iteration. A size1
correspond to a completely non-vectorized evaluation. A sizeNone
correspond to a completely vectorized one. - mpi_pool (
mpi_map.MPI_Pool
) – pool of processes to be used - import_set (set) – list of couples
(module_name,as_field)
to be imported asimport module_name as as_field
(for MPI purposes)
Returns: (
tuple
[\(N\)],list
)) – containing the \(N\) coefficients and log information from the optimizer.See also
TransportMaps.TriangularTransportMap.regression()
Note
the resulting coefficients \({\bf a}\) are automatically set at the end of the optimization. Use
coeffs()
in order to retrieve them.Note
The parameters
(qtype,qparams)
and(x,w)
are mutually exclusive, but one pair of them is necessary.- f (
-
reset_cache_minimize_kl_divergence_component
(cache)[source]¶ Reset the objective part of the cache space for the KL-divergence minimization
Parameters: params2 (dict) – dictionary of precomputed values
-
-
class
TransportMaps.Functionals.
MonotonicLinearSpanApproximation
(basis_list, spantype=None, order_list=None, multi_idxs=None, full_basis_list=None)[source]¶ Approximation of the type \(f \approx f_{\bf a} = \sum_{{\bf i} \in \mathcal{I}} {\bf a}_{\bf i} \Phi_{\bf i}\), monotonic in \(x_d\)
Parameters: -
minimize_kl_divergence_component
(f, x, w, x0=None, regularization=None, tol=0.0001, maxit=100, ders=2, fungrad=False, precomp_type='uni', batch_size=None, cache_level=1, mpi_pool=None)[source]¶ Compute \({\bf a}^\star = \arg\min_{\bf a}-\sum_{i=0}^m \log\pi\circ T_k(x_i) + \log\partial_{x_k}T_k(x_i) = \arg\min_{\bf a}-\sum_{i=0}^m f(x_i)\)
Parameters: - f (ProductDistributionParametricPullbackComponentFunction) – function \(f\)
- x (
ndarray
[\(m,d\)]) – quadrature points - w (
ndarray
[\(m\)]) – quadrature weights - x0 (
ndarray
[\(N\)]) – coefficients to be used as initial values for the optimization - regularization (dict) – defines the regularization to be used.
If
None
, no regularization is applied. If keytype=='L2'
then applies Tikonhov regularization with coefficient in keyalpha
. - tol (float) – tolerance to be used to solve the KL-divergence problem.
- maxit (int) – maximum number of iterations
- ders (int) – order of derivatives available for the solution of the optimization problem. 0 -> derivative free, 1 -> gradient, 2 -> hessian.
- fungrad (bool) – whether the distributions \(\pi_1,\pi_2\) provide the method
Distribution.tuple_grad_x_log_pdf()
computing the evaluation and the gradient in one step. This is used only forders==1
. - precomp_type (str) – whether to precompute univariate Vandermonde matrices ‘uni’ or multivariate Vandermonde matrices ‘multi’
- batch_size (
list
[3 or 2] ofint
orlist
ofbatch_size
) – the list contains the size of the batch to be used for each iteration. A size1
correspond to a completely non-vectorized evaluation. A sizeNone
correspond to a completely vectorized one. If the target distribution is aProductDistribution
, then the optimization problem decouples andbatch_size
is a list of lists containing the batch sizes to be used for each component of the map. - cache_level (int) – use high-level caching during the optimization, storing the
function evaluation
0
, and the gradient evaluation1
or nothing-1
- mpi_pool (
mpi_map.MPI_Pool
orlist
ofmpi_pool
) – pool of processes to be used,None
stands for one process. If the target distribution is aProductDistribution
, then the minimization problem decouples andmpi_pool
is a list containing ``mpi_pool``s for each component of the map.
-
precomp_regression
(x, precomp=None, *args, **kwargs)[source]¶ Precompute necessary structures for the speed up of
regression()
Parameters: Returns: (
dict
) – dictionary of necessary strucutres
-
regression
(f, fparams=None, d=None, qtype=None, qparams=None, x=None, w=None, x0=None, regularization=None, tol=0.0001, maxit=100, batch_size=(None, None), mpi_pool=None, import_set=set())[source]¶ Compute \({\bf a}^* = \arg\min_{\bf a} \Vert f - f_{\bf a} \Vert_{\pi}\).
Parameters: - f (
Function
orndarray
[\(m\)]) – function \(f\) or its functions values - d (Distribution) – distribution \(\pi\)
- fparams (dict) – parameters for function \(f\)
- qtype (int) – quadrature type to be used for the approximation of \(\mathbb{E}_{\pi}\)
- qparams (object) – parameters necessary for the construction of the quadrature
- x (
ndarray
[\(m,d\)]) – quadrature points used for the approximation of \(\mathbb{E}_{\pi}\) - w (
ndarray
[\(m\)]) – quadrature weights used for the approximation of \(\mathbb{E}_{\pi}\) - x0 (
ndarray
[\(N\)]) – coefficients to be used as initial values for the optimization - regularization (dict) – defines the regularization to be used.
If
None
, no regularization is applied. If keytype=='L2'
then applies Tikonhov regularization with coefficient in keyalpha
. - tol (float) – tolerance to be used to solve the regression problem.
- maxit (int) – maximum number of iterations
- batch_size (
list
[2] ofint
) – the list contains the size of the batch to be used for each iteration. A size1
correspond to a completely non-vectorized evaluation. A sizeNone
correspond to a completely vectorized one. - mpi_pool (
mpi_map.MPI_Pool
) – pool of processes to be used - import_set (set) – list of couples
(module_name,as_field)
to be imported asimport module_name as as_field
(for MPI purposes)
Returns: (
tuple
[\(N\)],list
)) – containing the \(N\) coefficients and log information from the optimizer.See also
TransportMaps.TriangularTransportMap.regression()
Note
the resulting coefficients \({\bf a}\) are automatically set at the end of the optimization. Use
coeffs()
in order to retrieve them.Note
The parameters
(qtype,qparams)
and(x,w)
are mutually exclusive, but one pair of them is necessary.- f (
-
-
class
TransportMaps.Functionals.
MonotonicIntegratedExponentialApproximation
(c, h, integ_ord_mult=6)[source]¶ Integrated Exponential approximation.
For \({\bf x} \in \mathbb{R}^d\) The approximation takes the form:
(1)¶\[f_{\bf a}({\bf x}) = c({\bf x};{\bf a}^c) + \int_0^{{\bf x}_d} \exp\left( h({\bf x}_{1:d-1},t;{\bf a}^e) \right) dt\]where
\[c({\bf x};{\bf a}^c) = \Phi({\bf x}) {\bf a}^c = \sum_{{\bf i}\in \mathcal{I}_c} \Phi_{\bf i}({\bf x}) {\bf a}^c_{\bf i} \qquad \text{and} \qquad h({\bf x}_{1:d-1},t;{\bf a}^e) = \Psi({\bf x}_{1:d-1},t) {\bf a}^e = \sum_{{\bf i}\in \mathcal{I}_e} \Psi_{\bf i}({\bf x}_{1:d-1},t) {\bf a}^e_{\bf i}\]for the set of basis \(\Phi\) and \(\Psi\) with cardinality \(\sharp \mathcal{I}_c = N_c\) and \(\sharp \mathcal{I}_e = N_e\). In the following \(N=N_c+N_e\).
Parameters: - c (
LinearSpanApproximation
) – \(d-1\) dimensional approximation of \(c({\bf x}_{1:d-1};{\bf a}^c)\). - h (
LinearSpanApproximation
) – \(d\) dimensional approximation of \(h({\bf x}_{1:d-1},t;{\bf a}^e)\). - integ_ord_mult (int) – multiplier for the number of Gauss points to be used
in the approximation of \(\int_0^{{\bf x}_d}\). The resulting number of
points is given by the product of the order in the \(d\) direction
and
integ_ord_mult
.
-
evaluate
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: (
ndarray
[\(m\)]) – function evaluations- x (
-
grad_a
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a} f_{\bf a}({\bf x})\)
- x (
-
grad_a_grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf a} \nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d\)]) – \(\nabla_{\bf a} \nabla_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_a_grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf a} \nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d\)]) – \(\nabla_{\bf a} \nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
grad_a_hess_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf a} \nabla^2_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d,d\)]) – \(\nabla_{\bf a} \nabla^2_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_a_hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf a}\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla_{\bf a}\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
grad_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
- x (
-
grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
hess_a
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla^2_{\bf a} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a} f_{\bf a}({\bf x})\)
- x (
-
hess_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cache
Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
hess_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d\)]) – \(\nabla^2_{\bf x} f_{\bf a}({\bf x})\)
- (
-
hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
partial2_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\partial^2_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m\)]) – \(\partial^2_{x_d} f_{\bf a}({\bf x})\)
- (
-
partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
precomp_Vandermonde_evaluate
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_Vandermonde_grad_x
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_Vandermonde_grad_x_partial_xd
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with the necessary structures
-
precomp_Vandermonde_hess_x
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\nabla^2_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_Vandermonde_hess_x_partial_xd
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with the necessary structures
-
precomp_Vandermonde_partial2_xd
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\partial^2_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with necessary structures
-
precomp_Vandermonde_partial_xd
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with necessary structures
-
precomp_evaluate
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_grad_x
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_grad_x_partial_xd
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with the necessary structures
-
precomp_hess_x
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\nabla^2_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_hess_x_partial_xd
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with the necessary structures
-
precomp_partial2_xd
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\partial^2_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with necessary structures
- c (
-
class
TransportMaps.Functionals.
MonotonicIntegratedSquaredApproximation
(c, h)[source]¶ Integrated Squared approximation.
For \({\bf x} \in \mathbb{R}^d\) The approximation takes the form:
\[f_{\bf a}({\bf x}) = c({\bf x};{\bf a}^c) + \int_0^{{\bf x}_d} \left( h({\bf x}_{1:d-1},t;{\bf a}^e) \right)^2 dt\]where
\[c({\bf x};{\bf a}^c) = \Phi({\bf x}) {\bf a}^c = \sum_{{\bf i}\in \mathcal{I}_c} \Phi_{\bf i}({\bf x}) {\bf a}^c_{\bf i} \qquad \text{and} \qquad h({\bf x}_{1:d-1},t;{\bf a}^e) = \Psi({\bf x}_{1:d-1},t) {\bf a}^e = \sum_{{\bf i}\in \mathcal{I}_e} \Psi_{\bf i}({\bf x}_{1:d-1},t) {\bf a}^e_{\bf i}\]for the set of basis \(\Phi\) and \(\Psi\) with cardinality \(\sharp \mathcal{I}_c = N_c\) and \(\sharp \mathcal{I}_e = N_e\). In the following \(N=N_c+N_e\).
Parameters: - c (
LinearSpanApproximation
) – \(d-1\) dimensional approximation of \(c({\bf x}_{1:d-1};{\bf a}^c)\). - h (
LinearSpanApproximation
) – \(d\) dimensional approximation of \(h({\bf x}_{1:d-1},t;{\bf a}^e)\).
-
evaluate
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: (
ndarray
[\(m\)]) – function evaluations- x (
-
grad_a
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a} f_{\bf a}({\bf x})\)
- x (
-
grad_a_grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla{\bf a} \nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d\)]) – \(\nabla{\bf a} \nabla_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_a_grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla{\bf a} \nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d\)]) – \(\nabla{\bf a} \nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
grad_a_hess_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla{\bf a} \nabla^2_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d,d\)]) – \(\nabla{\bf a} \nabla^2_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_a_hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla{\bf a} \nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,N,d,d\)]) – \(\nabla{\bf a} \nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
grad_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N\)]) – \(\nabla_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
grad_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
- (
-
grad_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d\)]) – \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
hess_a
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla^2_{\bf a} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a} f_{\bf a}({\bf x})\)
- x (
-
hess_a_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m,N,N\)]) – \(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
hess_x
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x} f_{\bf a}({\bf x})\)
- (
-
hess_x_partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m,d,d\)]) – \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)
- (
-
partial2_xd
(x, precomp=None, idxs_slice=slice(None, None, None), *args, **kwargs)[source]¶ Evaluate \(\partial^2_{x_d} f_{\bf a}\) at
x
.Parameters: Returns: - (
ndarray
[\(m\)]) – \(\partial^2_{x_d} f_{\bf a}({\bf x})\)
- (
-
partial_xd
(x, precomp=None, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\partial_{x_d} f_{\bf a}\) at
x
.Parameters: - x (
ndarray
[\(m,d\)]) – evaluation points - precomp (
dict
) – dictionary of precomputed values - idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (
dict
) – cache
Returns: - (
ndarray
[\(m\)]) – \(\partial_{x_d} f_{\bf a}({\bf x})\)
- x (
-
precomp_Vandermonde_evaluate
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_Vandermonde_grad_x
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_Vandermonde_grad_x_partial_xd
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with the necessary structures
-
precomp_Vandermonde_hess_x
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\nabla^2_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_Vandermonde_hess_x_partial_xd
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with the necessary structures
-
precomp_Vandermonde_partial2_xd
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\partial^2_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with necessary structures
-
precomp_Vandermonde_partial_xd
(x, precomp=None)[source]¶ Precompute necessary multi-variate structures for the evaluation of \(\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with necessary structures
-
precomp_evaluate
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_grad_x
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_grad_x_partial_xd
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with the necessary structures
-
precomp_hess_x
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\nabla^2_{\bf x} f_{\bf a}\) at
x
Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary containing the necessary structures
-
precomp_hess_x_partial_xd
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with the necessary structures
-
precomp_partial2_xd
(x, precomp=None, precomp_type='uni')[source]¶ Precompute necessary uni/multi-variate structures for the evaluation of \(\partial^2_{x_d} f_{\bf a}\) at
x
.Enriches the
precomp
dictionary if necessary.Parameters: Returns: (
dict
) – dictionary with necessary structures
- c (
-
class
TransportMaps.Functionals.
ProductDistributionParametricPullbackComponentFunction
(transport_map_component, base_distribution)[source]¶ Parametric function \(f[{\bf a}](x_{1:k}) = \log\pi\circ T_k[{\bf a}](x_{1:k}) + \log\partial_{x_k}T_k[{\bf a}](x_{1:k})\)
Parameters: - transport_map_component (MonotonicFunctionApproximation) – component \(T_k\) of monotone map \(T\)
- base_distribution (Distribution) – distribution \(\pi\)
-
coeffs
¶ Get the coefficients \({\bf a}\) of the function
-
evaluate
(x, params={}, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(f[{\bf a}](x_{1:k})\)
Parameters: - x (
ndarray
[\(m,k\)]) – evaluation points - params (dict) – parameters with keys
params_pi
,params_t
- idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cached values
Returns: (
ndarray
[\(m\)]) – evaluations- x (
-
grad_a
(x, params={}, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla_{\bf a}f[{\bf a}](x_{1:k})\)
Parameters: - x (
ndarray
[\(m,k\)]) – evaluation points - params (dict) – parameters with keys
params_pi
,params_t
- idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cached values
Returns: (
ndarray
[\(m,n\)]) – evaluations- x (
-
grad_a_hess_x
(x, params={}, idxs_slice=slice(None, None, None))[source]¶ Evaluate \(\nabla_{\bf a} \nabla^2_{\bf x} \log T_{k}^\sharp \pi({\bf x_{1:k}})\)
Parameters: Returns: - (
ndarray
[\(m,n,k,k\)]) – values of \(\nabla^2_{\bf x} \log T_{k}^\sharp \pi\) at the
x
points.
- (
-
grad_x
(x, params={}, idxs_slice=slice(None, None, None))[source]¶ Evaluate \(\nabla_{\bf x} \log T_{k}^\sharp \pi({\bf x_{1:k}})\)
Parameters: Returns: - (
ndarray
[\(m\)]) – values of \(\nabla_{\bf x} \log T_{k}^\sharp \pi\) at the
x
points.
- (
-
hess_a
(x, params={}, idxs_slice=slice(None, None, None), cache=None)[source]¶ Evaluate \(\nabla^2_{\bf a}f[{\bf a}](x_{1:k})\)
Parameters: - x (
ndarray
[\(m,k\)]) – evaluation points - params (dict) – parameters with keys
params_pi
,params_t
- idxs_slice (slice) – if precomputed values are present, this parameter
indicates at which of the points to evaluate. The number of indices
represented by
idxs_slice
must matchx.shape[0]
. - cache (dict) – cached values
Returns: (
ndarray
[\(m,n\)]) – evaluations- x (
-
hess_x
(x, params={}, idxs_slice=slice(None, None, None))[source]¶ Evaluate \(\nabla^2_{\bf x} \log T_{k}^\sharp \pi({\bf x_{1:k}})\)
Parameters: Returns: - (
ndarray
[\(m,k,k\)]) – values of \(\nabla^2_{\bf x} \log T_{k}^\sharp \pi\) at the
x
points.
- (
-
n_coeffs
¶ Get the number \(N\) of coefficients
-
class
TransportMaps.Functionals.
MonotonicFrozenFunction
(dim)[source]¶ [Abstract] Frozen function. No optimization over the coefficients allowed.
-
precomp_evaluate
(x, *args, **kwargs)[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_grad_x
(x, *args, **kwargs)[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\nabla_{\bf x} f_{\bf a}\) at
x
Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
precomp_grad_x_partial_xd
(x, *args, **kwargs)[source]¶ [Abstract] Precompute necessary structures for the evaluation of \(\nabla_{\bf x}\partial_{x_d} f_{\bf a}\) at
x
.Parameters: x ( ndarray
[\(m,d\)]) – evaluation pointsReturns: ( dict
) – data structures
-
-
class
TransportMaps.Functionals.
FrozenLinear
(dim, a1, a2)[source]¶ Frozen Linear map \({\bf x} \rightarrow a_1 + a_2 {\bf x}_d\)
Parameters: -
evaluate
(x, *args, **kwargs)[source]¶ Evaluate \(f_{\bf a}({\bf x}) = a_1 + a_2 {\bf x}_d\)
Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: ( ndarray
[\(m\)]) – function value at points x
-
grad_x
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
This is:
\[\begin{split}\nabla_{\bf x} f_{\bf a}({\bf x}) = \begin{bmatrix} \partial_{{\bf x}_1} f_{\bf a}({\bf x}) \\ \partial_{{\bf x}_2} f_{\bf a}({\bf x}) \\ \vdots \\ \partial_{{\bf x}_d} f_{\bf a}({\bf x}) \end{bmatrix} = \begin{bmatrix} 0 \\ 0 \\ \vdots \\ a2 \end{bmatrix}\end{split}\]Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d\)]) – - \(\nabla_{\bf x}f_{\bf a}({\bf x};{\bf a})\) at points x
- (
-
grad_x_partial_xd
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x}\partial_{{\bf x}_d} f_{\bf a}({\bf x}) = 0\)
Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d\)]) – - \(\nabla_{\bf x}\partial_{{\bf x}_d}f_{\bf a}({\bf x})\) at points x
- (
-
inverse
(xkm1, y)[source]¶ Compute \(f_{\bf a}^{-1}({\bf x}_{1:d-1})(y):={\bf x}_d\) s.t. \(f_{\bf a}({\bf x}_{1:d-1},{\bf x}_d) - y = 0\).
Due to the form of the approximation we have:
\[f_{\bf a}^{-1}({\bf x}_{1:d-1})(y) = \frac{y-a_1}{a_2}\]Parameters: Returns: (
float
) – inverse value \(x\).
-
n_coeffs
¶ Returns – 2
-
partial2_xd
(x, *args, **kwargs)[source]¶ Evaluate \(\partial^2_{{\bf x}_d} f_{\bf a}({\bf x}) = 0\)
Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m\)]) – - \(\partial^2_{{\bf x}_d}f_{\bf a}({\bf x};{\bf a})\) at points x
- (
-
-
class
TransportMaps.Functionals.
FrozenExponential
(dim)[source]¶ Frozen Exponential map \(f_{\bf a}:{\bf x} \mapsto \exp( {\bf x}_d )\)
Parameters: dim (int) – input dimension \(d\) -
evaluate
(x, *args, **kwargs)[source]¶ Evaluate \(f_{\bf a}({\bf x}) = \exp({\bf x}_d)\)
Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: ( ndarray
[\(m\)]) – function value at points x
-
grad_x
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
This is:
\[\begin{split}\nabla_{\bf x} f_{\bf a}({\bf x}) = \begin{bmatrix} \partial_{{\bf x}_1} f_{\bf a}({\bf x}) \\ \partial_{{\bf x}_2} f_{\bf a}({\bf x}) \\ \vdots \\ \partial_{{\bf x}_d} f_{\bf a}({\bf x}) \end{bmatrix} = \begin{bmatrix} 0 \\ 0 \\ \vdots \\ \exp({\bf x}_d) \end{bmatrix}\end{split}\]Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d\)]) – - \(\nabla_{\bf x}f_{\bf a}({\bf x};{\bf a})\) at points x
- (
-
grad_x_partial_xd
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x}\partial_{{\bf x}_d} f_{\bf a}({\bf x})\)
This is:
\[\begin{split}\nabla_{\bf x} \partial_{{\bf x}_d} f_{\bf a}({\bf x}) = \begin{bmatrix} \partial_{{\bf x}_1} f_{\bf a}({\bf x}) \\ \partial_{{\bf x}_2} f_{\bf a}({\bf x}) \\ \vdots \\ \partial_{{\bf x}_d} f_{\bf a}({\bf x}) \end{bmatrix} = \begin{bmatrix} 0 \\ 0 \\ \vdots \\ \exp({\bf x}_d) \end{bmatrix}\end{split}\]Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d\)]) – - \(\nabla_{\bf x}\partial_{{\bf x}_d}f_{\bf a}({\bf x})\) at points x
- (
-
hess_x
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x} f_{\bf a}({\bf x})\)
This is:
\[\begin{split}\nabla^2_{\bf x} f_{\bf a}({\bf x}) = \begin{bmatrix} \partial^2_{{\bf x}_1} f_{\bf a}({\bf x}) & \partial_{{\bf x}_1{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial_{{\bf x}_1{\bf x}_d} f_{\bf a}({\bf x}) \\ \partial_{{\bf x}_2 {\bf x}_1} f_{\bf a}({\bf x}) & \partial^2_{{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial_{{\bf x}_2{\bf x}_d} f_{\bf a}({\bf x}) \\ \vdots & & \ddots & \\ \partial_{{\bf x}_d{\bf x}_1} f_{\bf a}({\bf x}) & \partial_{{\bf x}_d{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial^2_{{\bf x}_d} f_{\bf a}({\bf x}) \end{bmatrix} = \begin{bmatrix} 0 & \cdots & 0 & 0 \\ \vdots & \ddots & 0 & \vdots \\ 0 & \cdots & 0 & 0 \\ 0 & \cdots & 0 & \exp({\bf x}_d) \end{bmatrix}\end{split}\]Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d,d\)]) – - \(\nabla^2_{\bf x}f_{\bf a}({\bf x})\) at points x
- (
-
inverse
(xkm1, y)[source]¶ Compute \(f_{\bf a}^{-1}({\bf x}_{1:d-1})(y):={\bf x}_d\) s.t. \(f_{\bf a}({\bf x}_{1:d-1},{\bf x}_d) - y = 0\).
Due to the form of the approximation we have:
\[f_{\bf a}^{-1}({\bf x}_{1:d-1})(y) = \log(y)\]Parameters: Returns: (
float
) – inverse value \(x\).
-
n_coeffs
¶ Returns – 0
-
partial2_xd
(x, *args, **kwargs)[source]¶ Evaluate \(\partial^2_{{\bf x}_d} f_{\bf a}({\bf x}) = \exp({\bf x}_d)\)
Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m\)]) – - \(\partial^2_{{\bf x}_d}f_{\bf a}({\bf x})\) at points x
- (
-
-
class
TransportMaps.Functionals.
FrozenGaussianToUniform
(dim)[source]¶ Frozen Gaussian To Uniform map.
This is given by the Cumulative Distribution Function of a standard normal distribution along the last coordinate:
\[f_{\bf a}({\bf x}) = \frac{1}{2} \left[ 1 + \text{erf}\left( \frac{x}{\sqrt{2}} \right)\right]\]Parameters: dim (int) – input dimension \(d\) -
evaluate
(x, *args, **kwargs)[source]¶ Evaluate \(f_{\bf a}({\bf x})\)
Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: ( ndarray
[\(m\)]) – function value at points x
-
grad_x
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x} f_{\bf a}({\bf x})\)
This is:
\[\begin{split}\nabla_{\bf x} f_{\bf a}({\bf x}) = \begin{bmatrix} \partial_{{\bf x}_1} f_{\bf a}({\bf x}) \\ \partial_{{\bf x}_2} f_{\bf a}({\bf x}) \\ \vdots \\ \partial_{{\bf x}_d} f_{\bf a}({\bf x}) \end{bmatrix} = \begin{bmatrix} 0 \\ 0 \\ \vdots \\ (2\pi)^{-1}\exp(-\frac{{\bf x}^2_d}{2}) \end{bmatrix}\end{split}\]Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d\)]) – - \(\nabla_{\bf x}f_{\bf a}({\bf x})\) at points x
- (
-
grad_x_partial_xd
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla_{\bf x}\partial_{{\bf x}_d} f_{\bf a}({\bf x})\)
This is:
\[\begin{split}\nabla_{\bf x} \partial_{{\bf x}_d} f_{\bf a}({\bf x}) = \begin{bmatrix} \partial_{{\bf x}_1} f_{\bf a}({\bf x}) \\ \partial_{{\bf x}_2} f_{\bf a}({\bf x}) \\ \vdots \\ \partial_{{\bf x}_d} f_{\bf a}({\bf x}) \end{bmatrix} = \begin{bmatrix} 0 \\ 0 \\ \vdots \\ -{\bf x}_d (2\pi)^{-1} \exp(-\frac{{\bf x}^2_d}{2}) \end{bmatrix}\end{split}\]Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d\)]) – - \(\nabla_{\bf x}\partial_{{\bf x}_d}f_{\bf a}({\bf x})\) at points x
- (
-
hess_x
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x} f_{\bf a}({\bf x})\)
This is:
\[\begin{split}\nabla^2_{\bf x} f_{\bf a}({\bf x}) = \begin{bmatrix} \partial^2_{{\bf x}_1} f_{\bf a}({\bf x}) & \partial_{{\bf x}_1{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial_{{\bf x}_1{\bf x}_d} f_{\bf a}({\bf x}) \\ \partial_{{\bf x}_2 {\bf x}_1} f_{\bf a}({\bf x}) & \partial^2_{{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial_{{\bf x}_2{\bf x}_d} f_{\bf a}({\bf x}) \\ \vdots & & \ddots & \\ \partial_{{\bf x}_d{\bf x}_1} f_{\bf a}({\bf x}) & \partial_{{\bf x}_d{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial^2_{{\bf x}_d} f_{\bf a}({\bf x}) \end{bmatrix} = \begin{bmatrix} 0 & \cdots & 0 & 0 \\ \vdots & \ddots & 0 & \vdots \\ 0 & \cdots & 0 & 0 \\ 0 & \cdots & 0 & -{\bf x}_d (2\pi)^{-1} \exp(-\frac{{\bf x}^2_d}{2}) \end{bmatrix}\end{split}\]Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d,d\)]) – - \(\nabla^2_{\bf x}f_{\bf a}({\bf x})\) at points x
- (
-
hess_x_partial_xd
(x, *args, **kwargs)[source]¶ Evaluate \(\nabla^2_{\bf x} \partial_{{\bf x}_d} f_{\bf a}({\bf x})\)
This is:
\[\begin{split}\nabla^2_{\bf x} \partial_{{\bf x}_d} f_{\bf a}({\bf x}) = \begin{bmatrix} \partial^2_{{\bf x}_1} f_{\bf a}({\bf x}) & \partial_{{\bf x}_1{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial_{{\bf x}_1{\bf x}_d} f_{\bf a}({\bf x}) \\ \partial_{{\bf x}_2 {\bf x}_1} f_{\bf a}({\bf x}) & \partial^2_{{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial_{{\bf x}_2{\bf x}_d} f_{\bf a}({\bf x}) \\ \vdots & & \ddots & \\ \partial_{{\bf x}_d{\bf x}_1} f_{\bf a}({\bf x}) & \partial_{{\bf x}_d{\bf x}_2} f_{\bf a}({\bf x}) & \cdots & \partial^2_{{\bf x}_d} f_{\bf a}({\bf x}) \end{bmatrix} = \begin{bmatrix} 0 & \cdots & 0 & 0 \\ \vdots & \ddots & 0 & \vdots \\ 0 & \cdots & 0 & 0 \\ 0 & \cdots & 0 & ({\bf x}_d - 1) (2\pi)^{-1} \exp(-\frac{{\bf x}^2_d}{2}) \end{bmatrix}\end{split}\]Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m,d,d\)]) – - \(\nabla^2_{\bf x}\partial_{{\bf x}_d}f_{\bf a}({\bf x})\) at points x
- (
-
invserse
(xkm1, y)[source]¶ Compute \(f_{\bf a}^{-1}({\bf x}_{1:d-1})(y):={\bf x}_d\) s.t. \(f_{\bf a}({\bf x}_{1:d-1},{\bf x}_d) - y = 0\).
Due to the form of the approximation we have:
\[f_{\bf a}^{-1}({\bf x}_{1:d-1})(y) = \sqrt{2} \text{erf}^{-1}(2y-1)\]Parameters: Returns: (
float
) – inverse value \(x\).
-
n_coeffs
¶ Returns – 0
-
partial2_xd
(x, *args, **kwargs)[source]¶ Evaluate \(\partial^2_{{\bf x}_d} f_{\bf a}({\bf x}) = -{\bf x}_d (2\pi)^{-1} \exp(-\frac{{\bf x}^2_d}{2})\)
Parameters: x ( ndarray
[\(m,d\)]) – evaluation points.Returns: - (
ndarray
[\(m\)]) – - \(\partial^2_{{\bf x}_d}f_{\bf a}({\bf x})\) at points x
- (
-