TransportMaps.Maps.Functionals.IntegratedSquaredParametricMonotoneFunctionalBase

Module Contents

Classes

IntegratedSquaredParametricMonotoneFunctional

Integrated Squared approximation.

MonotonicIntegratedSquaredApproximation

Integrated Squared approximation.

class TransportMaps.Maps.Functionals.IntegratedSquaredParametricMonotoneFunctionalBase.IntegratedSquaredParametricMonotoneFunctional(c, h)[source]

Bases: TransportMaps.Maps.Functionals.ParametricMonotoneFunctionalBase.ParametricMonotoneFunctional

Integrated Squared 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} \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 (LinearSpanTensorizedParametricFunctional) – \(d-1\) dimensional approximation of \(c({\bf x}_{1:d-1};{\bf a}^c)\).

  • h (LinearSpanTensorizedParametricFunctional) – \(d\) dimensional approximation of \(h({\bf x}_{1:d-1},t;{\bf a}^e)\).

property n_coeffs[source]

Get the number \(N\) of coefficients \({\bf a}\)

Returns:

(int) – number of coefficients

property coeffs[source]

Get the coefficients \({\bf a}\)

Returns:

(ndarray [\(N\)]) – coefficients

get_ncalls_tree(indent='')[source]
get_nevals_tree(indent='')[source]
get_teval_tree(indent='')[source]
update_ncalls_tree(obj)[source]
update_nevals_tree(obj)[source]
update_teval_tree(obj)[source]
reset_counters()[source]
init_coeffs()[source]

Initialize the coefficients \({\bf a}\)

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

  • precomp_type (str) – whether to precompute uni-variate Vandermonde matrices (uni) or to precompute the multi-variate Vandermonde matrices (multi)

Returns:

(dict) – dictionary containing the necessary structures

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

Returns:

(dict) – dictionary containing the necessary structures

evaluate(x, precomp=None, idxs_slice=slice(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 match x.shape[0].

  • cache (dict) – cache

Returns:

(ndarray [\(m,1\)]) – function evaluations

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

  • precomp_type (str) – whether to precompute uni-variate Vandermonde matrices (uni) or to precompute the multi-variate Vandermonde matrices (multi)

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

Returns:

(dict) – dictionary containing the necessary structures

grad_x(x, precomp=None, idxs_slice=slice(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 match x.shape[0].

Returns:

(ndarray [\(m,1,d\)]) –

\(\nabla_{\bf x} f_{\bf a}({\bf x})\)

grad_a_grad_x(x, precomp=None, idxs_slice=slice(None), *args, **kwargs)[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 match x.shape[0].

Returns:

(ndarray [\(m,1,N,d\)]) –

\(\nabla{\bf a} \nabla_{\bf x} f_{\bf a}({\bf x})\)

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

  • precomp_type (str) – whether to precompute uni-variate Vandermonde matrices (uni) or to precompute the multi-variate Vandermonde matrices (multi)

Returns:

(dict) – dictionary containing 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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

Returns:

(dict) – dictionary containing the necessary structures

hess_x(x, precomp=None, idxs_slice=slice(None), *args, **kwargs)[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 match x.shape[0].

Returns:

(ndarray [\(m,1,d,d\)]) –

\(\nabla^2_{\bf x} f_{\bf a}({\bf x})\)

grad_a_hess_x(x, precomp=None, idxs_slice=slice(None), *args, **kwargs)[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 match x.shape[0].

Returns:

(ndarray [\(m,1,N,d,d\)]) –

\(\nabla{\bf a} \nabla^2_{\bf x} f_{\bf a}({\bf x})\)

grad_a(x, precomp=None, idxs_slice=slice(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 match x.shape[0].

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,N\)]) –

\(\nabla_{\bf a} f_{\bf a}({\bf x})\)

hess_a(x, precomp=None, idxs_slice=slice(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 match x.shape[0].

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,N,N\)]) –

\(\nabla^2_{\bf a} f_{\bf a}({\bf x})\)

precomp_partial_xd(x, precomp=None, precomp_type='uni')[source]

Precompute necessary uni/multi-variate structures for the evaluation of \(\partial_{x_d} f_{\bf a}\) at x.

Enriches the precomp dictionary if necessary.

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

  • precomp_type (str) – whether to precompute uni-variate Vandermonde matrices (uni) or to precompute the multi-variate Vandermonde matrices (multi)

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

Returns:

(dict) – dictionary with necessary structures

partial_xd(x, precomp=None, idxs_slice=slice(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 match x.shape[0].

  • cache (dict) – cache

Returns:

(ndarray [\(m,1\)]) –

\(\partial_{x_d} f_{\bf a}({\bf x})\)

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

  • precomp_type (str) – whether to precompute uni-variate Vandermonde matrices (uni) or to precompute the multi-variate Vandermonde matrices (multi)

Returns:

(dict) – dictionary with 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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

Returns:

(dict) – dictionary with the necessary structures

grad_x_partial_xd(x, precomp=None, idxs_slice=slice(None), *args, **kwargs)[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 match x.shape[0].

Returns:

(ndarray [\(m,1,d\)]) –

\(\nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)

grad_a_grad_x_partial_xd(x, precomp=None, idxs_slice=slice(None), *args, **kwargs)[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 match x.shape[0].

Returns:

(ndarray [\(m,1,N,d\)]) –

\(\nabla{\bf a} \nabla_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

  • precomp_type (str) – whether to precompute uni-variate Vandermonde matrices (uni) or to precompute the multi-variate Vandermonde matrices (multi)

Returns:

(dict) – dictionary with 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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

Returns:

(dict) – dictionary with the necessary structures

hess_x_partial_xd(x, precomp=None, idxs_slice=slice(None), *args, **kwargs)[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 match x.shape[0].

Returns:

(ndarray [\(m,1,d,d\)]) –

\(\nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)

grad_a_hess_x_partial_xd(x, precomp=None, idxs_slice=slice(None), *args, **kwargs)[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 match x.shape[0].

Returns:

(ndarray [\(m,1,N,d,d\)]) –

\(\nabla{\bf a} \nabla^2_{\bf x}\partial_{x_d} f_{\bf a}({\bf x})\)

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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

  • precomp_type (str) – whether to precompute uni-variate Vandermonde matrices (uni) or to precompute the multi-variate Vandermonde matrices (multi)

Returns:

(dict) – dictionary with 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:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary of precomputed values

Returns:

(dict) – dictionary with necessary structures

partial2_xd(x, precomp=None, idxs_slice=slice(None), *args, **kwargs)[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 match x.shape[0].

Returns:

(ndarray [\(m,1\)]) –

\(\partial^2_{x_d} f_{\bf a}({\bf x})\)

grad_a_partial_xd(x, precomp=None, idxs_slice=slice(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 match x.shape[0].

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,N\)]) –

\(\nabla_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)

hess_a_partial_xd(x, precomp=None, idxs_slice=slice(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 match x.shape[0].

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,N,N\)]) –

\(\nabla^2_{\bf a}\partial_{x_d} f_{\bf a}({\bf x})\)

get_identity_coeffs()[source]
precomp_regression(x, precomp=None, *args, **kwargs)[source]

Precompute necessary structures for the speed up of regression()

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • precomp (dict) – dictionary to be updated

Returns:

(dict) – dictionary of necessary strucutres

class TransportMaps.Maps.Functionals.IntegratedSquaredParametricMonotoneFunctionalBase.MonotonicIntegratedSquaredApproximation(*args, **kwars)[source]

Bases: IntegratedSquaredParametricMonotoneFunctional

Integrated Squared approximation.

For \({\bf x} \in \mathbb{R}^d\) The approximation takes the form:

(2)\[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 (LinearSpanTensorizedParametricFunctional) – \(d-1\) dimensional approximation of \(c({\bf x}_{1:d-1};{\bf a}^c)\).

  • h (LinearSpanTensorizedParametricFunctional) – \(d\) dimensional approximation of \(h({\bf x}_{1:d-1},t;{\bf a}^e)\).