TransportMaps.L2.minimize_L2_distance

Module Contents

Functions

map_regression(tm, t[, tparams, d, qtype, qparams, x, ...])

Compute \({\bf a}^* = \arg\min_{\bf a} \Vert T - T[{\bf a}] \Vert_{\pi}\).

functional_regression(fn, f[, fparams, d, qtype, ...])

Compute \({\bf a}^* = \arg\min_{\bf a} \Vert f - f_{\bf a} \Vert_{\pi}\).

functional_regression_objective(a, params)

Objective function \(\Vert f - f_{\bf a} \Vert^2_{\pi}\)

functional_regression_grad_a_objective(a, params)

Objective function \(\nabla_{\bf a} \Vert f - f_{\bf a} \Vert^2_{\pi}\)

functional_regression_hess_a_objective(a, params)

Objective function \(\nabla_{\bf a}^2 \Vert f - f_{\bf a} \Vert^2_{\pi}\)

functional_regression_action_storage_hess_a_objective(a, ...)

Assemble/fetch Hessian \(\nabla_{\bf a}^2 \Vert f - f_{\bf a} \Vert^2_{\pi}\) and evaluate its action on \(v\)

affine_triangular_map_regression(tm, t[, tparams, d, ...])

Compute \({\bf a}^* = \arg\min_{\bf a} \Vert T - T({\bf a}) \Vert_{\pi}\).

TransportMaps.L2.minimize_L2_distance.map_regression(tm: TransportMaps.Maps.ParametricComponentwiseMap, t, tparams=None, d=None, qtype=None, qparams=None, x=None, w=None, x0=None, regularization=None, tol=0.0001, maxit=100, batch_size_list=None, mpi_pool_list=None)[source]

Compute \({\bf a}^* = \arg\min_{\bf a} \Vert T - T[{\bf a}] \Vert_{\pi}\).

This regression problem can be completely decoupled if the measure is a product measure, obtaining

\[a^{(i)*} = \arg\min_{\bf a^{(i)}} \Vert T_i - T_i[{\bf a}^{(i)}] \Vert_{\pi_i}\]
Parameters:
  • tm (ParametricComponentwiseMap) – transport map \(T\)

  • t (function or ndarray [\(m\)]) – function \(t\) with signature t(x) or its functions values

  • tparams (dict) – parameters for function \(t\)

  • 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 key type=='L2' then applies Tikonhov regularization with coefficient in key alpha.

  • tol (float) – tolerance to be used to solve the regression problem.

  • maxit (int) – maximum number of iterations

  • batch_size_list (list [d] tuple [3] int) – Each of the tuples in the list corresponds to each component of the map. The entries of the tuple define whether to evaluate the regression in batches of a certain size or not. A size 1 correspond to a completely non-vectorized evaluation. A size None correspond to a completely vectorized one. (Note: if nprocs > 1, then the batch size defines the size of the batch for each process)

  • mpi_pool_list (list [d] mpi_map.MPI_Pool or None) – pool of processes to be used for function evaluation, gradient evaluation and Hessian evaluation for each component of the approximation. Value None will use serial evaluation.

Returns:

(list [\(d\)]) containing log information from the optimizer.

See also

MonotonicApproximation

Note

the resulting coefficients \({\bf a}\) are automatically set at the end of the optimization. Use get_coeffs() in order to retrieve them.

Note

The parameters (qtype,qparams) and (x,w) are mutually exclusive, but one pair of them is necessary.

TransportMaps.L2.minimize_L2_distance.functional_regression(fn: TransportMaps.Maps.Functionals.ParametricFunctional, 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:
  • fn (ParametricFunctional) – the function \(f_{\bf a}\)

  • f (Function or ndarray [\(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 key type=='L2' then applies Tikonhov regularization with coefficient in key alpha.

  • tol (float) – tolerance to be used to solve the regression problem.

  • maxit (int) – maximum number of iterations

  • batch_size (list [3] of int) – the list contains the size of the batch to be used for each iteration. A size 1 correspond to a completely non-vectorized evaluation. A size None 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 as import 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.

TransportMaps.L2.minimize_L2_distance.functional_regression_objective(a, params)[source]

Objective function \(\Vert f - f_{\bf a} \Vert^2_{\pi}\)

Parameters:
  • a (ndarray [\(N\)]) – coefficients

  • params (dict) – dictionary of parameters

TransportMaps.L2.minimize_L2_distance.functional_regression_grad_a_objective(a, params)[source]

Objective function \(\nabla_{\bf a} \Vert f - f_{\bf a} \Vert^2_{\pi}\)

Parameters:
  • a (ndarray [\(N\)]) – coefficients

  • params (dict) – dictionary of parameters

TransportMaps.L2.minimize_L2_distance.functional_regression_hess_a_objective(a, params)[source]

Objective function \(\nabla_{\bf a}^2 \Vert f - f_{\bf a} \Vert^2_{\pi}\)

Parameters:
  • a (ndarray [\(N\)]) – coefficients

  • params (dict) – dictionary of parameters

TransportMaps.L2.minimize_L2_distance.functional_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:
  • a (ndarray [\(N\)]) – coefficients

  • v (ndarray [\(N\)]) – vector on which to apply the Hessian

  • params (dict) – dictionary of parameters

TransportMaps.L2.minimize_L2_distance.affine_triangular_map_regression(tm: TransportMaps.Maps.ParametricComponentwiseMap, t, tparams=None, d=None, qtype=None, qparams=None, x=None, w=None, regularization=None, **kwargs)[source]

Compute \({\bf a}^* = \arg\min_{\bf a} \Vert T - T({\bf a}) \Vert_{\pi}\).

This regression problem can be completely decoupled if the measure is a product measure, obtaining

\[a^{(i)*} = \arg\min_{\bf a^{(i)}} \Vert T_i - T_i({\bf a}^{(i)}) \Vert_{\pi_i}\]
Parameters:
  • tm (ParametricComponentwiseMap) – map \(T({\bf a})\)

  • t (function or ndarray [\(m\)]) – function \(t\) with signature t(x) or its functions values

  • tparams (dict) – parameters for function \(t\)

  • 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}\)

  • regularization (dict) – defines the regularization to be used. If None, no regularization is applied. If key type=='L2' then applies Tikonhov regularization with coefficient in key alpha.

Returns:

(tuple (ndarray [\(N\)], list)) – containing the \(N\) coefficients and log information.

Note

the resulting coefficients \({\bf a}\) are automatically set at the end of the optimization. Use get_coeffs() in order to retrieve them.

Note

The parameters (qtype,qparams) and (x,w) are mutually exclusive, but one pair of them is necessary.