Functions

The pyunlocbox.functions module implements an interface for solvers to access the functions to be optimized as well as common objective functions.

Interface

The func base class defines a common interface to all functions:

func.cap(x)

Test the capabilities of the function object.

func.eval(x)

Function evaluation.

func.prox(x, T)

Function proximal operator.

func.grad(x)

Function gradient.

Functions

Then, derived classes implement various common objective functions.

Norm operators (based on norm)

norm_l1(**kwargs)

L1-norm (eval, prox).

norm_l2(**kwargs)

L2-norm (eval, prox, grad).

norm_nuclear(**kwargs)

Nuclear-norm (eval, prox).

norm_tv([dim, verbosity])

TV-norm (eval, prox).

Projection operators (based on proj)

proj_positive(**kwargs)

Projection on the positive octant (eval, prox).

proj_b2([epsilon, method])

Projection on the L2-ball (eval, prox).

proj_lineq([A, pinvA])

Projection on the plane satisfying the linear equality Az = y (eval, prox).

proj_spsd(**kwargs)

Projection on symmetric positive semi-definite matrices (eval, prox).

Miscellaneous

dummy(**kwargs)

Dummy function (eval, prox, grad).

structured_sparsity([lambda_, groups, weights])

Structured sparsity (eval, prox).

Inheritance diagram of pyunlocbox.functions
class pyunlocbox.functions.func(y=0, A=None, At=None, tight=True, nu=1, tol=0.001, maxit=200, **kwargs)[source]

Bases: object

This class defines the function object interface.

It is intended to be a base class for standard functions which will implement the required methods. It can also be instantiated by user code and dynamically modified for rapid testing. The instanced objects are meant to be passed to the pyunlocbox.solvers.solve() solving function.

Parameters
yarray_like, optional

Measurements. Default is 0.

Afunction or ndarray, optional

The forward operator. Default is the identity, \(A(x)=x\). If A is an ndarray, it will be converted to the operator form.

Atfunction or ndarray, optional

The adjoint operator. If At is an ndarray, it will be converted to the operator form. If A is an ndarray, default is the transpose of A. If A is a function, default is A, \(At(x)=A(x)\).

tightbool, optional

True if A is a tight frame (semi-orthogonal linear transform), False otherwise. Default is True.

nufloat, optional

Bound on the norm of the operator A, i.e. \(\|A(x)\|^2 \leq \nu \|x\|^2\). Default is 1.

tolfloat, optional

The tolerance stopping criterion. The exact definition depends on the function object, please see the documentation of the considered function. Default is 1e-3.

maxitint, optional

The maximum number of iterations. Default is 200.

Examples

Let’s define a parabola as an example of the manual implementation of a function object :

>>> from pyunlocbox import functions
>>> f = functions.func()
>>> f._eval = lambda x: x**2
>>> f._grad = lambda x: 2*x
>>> x = [1, 2, 3, 4]
>>> f.eval(x)
array([ 1,  4,  9, 16])
>>> f.grad(x)
array([2, 4, 6, 8])
>>> f.cap(x)
['EVAL', 'GRAD']
eval(x)[source]

Function evaluation.

Parameters
xarray_like

The evaluation point. If x is a matrix, the function gets evaluated for each column, as if it was a set of independent problems. Some functions, like the nuclear norm, are only defined on matrices.

Returns
zfloat

The objective function evaluated at x. If x is a matrix, the sum of the objectives is returned.

Notes

This method is required by the pyunlocbox.solvers.solve() solving function to evaluate the objective function. Each function class should therefore define it.

prox(x, T)[source]

Function proximal operator.

Parameters
xarray_like

The evaluation point. If x is a matrix, the function gets evaluated for each column, as if it was a set of independent problems. Some functions, like the nuclear norm, are only defined on matrices.

Tfloat

The regularization parameter.

Returns
zndarray

The proximal operator evaluated for each column of x.

Notes

The proximal operator is defined by \(\operatorname{prox}_{\gamma f}(x) = \operatorname{arg\,min} \limits_z \frac{1}{2} \|x-z\|_2^2 + \gamma f(z)\)

This method is required by some solvers.

When the map A in the function construction is a tight frame (semi-orthogonal linear transformation), we can use property (x) of Table 10.1 in [CP11] to compute the proximal operator of the composition of A with the base function. Whenever this is not the case, we have to resort to some iterative procedure, which may be very inefficient.

grad(x)[source]

Function gradient.

Parameters
xarray_like

The evaluation point. If x is a matrix, the function gets evaluated for each column, as if it was a set of independent problems. Some functions, like the nuclear norm, are only defined on matrices.

Returns
zndarray

The objective function gradient evaluated for each column of x.

Notes

This method is required by some solvers.

cap(x)[source]

Test the capabilities of the function object.

Parameters
xarray_like

The evaluation point. Not really needed, but this function calls the methods of the object to test if they can properly execute without raising an exception. Therefore it needs some evaluation point with a consistent size.

Returns
caplist of string

A list of capabilities (‘EVAL’, ‘GRAD’, ‘PROX’).

class pyunlocbox.functions.dummy(**kwargs)[source]

Bases: pyunlocbox.functions.func

Dummy function (eval, prox, grad).

This can be used as a second function object when there is only one function to minimize. It always evaluates as 0.

Examples

>>> from pyunlocbox import functions
>>> f = functions.dummy()
>>> x = [1, 2, 3, 4]
>>> f.eval(x)
0
>>> f.prox(x, 1)
array([1, 2, 3, 4])
>>> f.grad(x)
array([0, 0, 0, 0])
class pyunlocbox.functions.norm(lambda_=1, w=1, **kwargs)[source]

Bases: pyunlocbox.functions.func

Base class which defines the attributes of the norm objects.

See generic attributes descriptions of the pyunlocbox.functions.func base class.

Parameters
lambda_float, optional

Regularization parameter \(\lambda\). Default is 1.

warray_like, optional

Weights for a weighted norm. Default is 1.

class pyunlocbox.functions.norm_l1(**kwargs)[source]

Bases: pyunlocbox.functions.norm

L1-norm (eval, prox).

See generic attributes descriptions of the pyunlocbox.functions.norm base class. Note that the constructor takes keyword-only parameters.

Notes

  • The L1-norm of the vector x is given by \(\lambda \|w \cdot (A(x)-y)\|_1\).

  • The L1-norm proximal operator evaluated at x is given by \(\operatorname{arg\,min}\limits_z \frac{1}{2} \|x-z\|_2^2 + \gamma \|w \cdot (A(z)-y)\|_1\) where \(\gamma = \lambda \cdot T\). This is simply a soft thresholding.

Examples

>>> from pyunlocbox import functions
>>> f = functions.norm_l1()
>>> f.eval([1, 2, 3, 4])
10
>>> f.prox([1, 2, 3, 4], 1)
array([0., 1., 2., 3.])
class pyunlocbox.functions.norm_l2(**kwargs)[source]

Bases: pyunlocbox.functions.norm

L2-norm (eval, prox, grad).

See generic attributes descriptions of the pyunlocbox.functions.norm base class. Note that the constructor takes keyword-only parameters.

Notes

  • The squared L2-norm of the vector x is given by \(\lambda \|w \cdot (A(x)-y)\|_2^2\).

  • The squared L2-norm proximal operator evaluated at x is given by \(\operatorname{arg\,min}\limits_z \frac{1}{2} \|x-z\|_2^2 + \gamma \|w \cdot (A(z)-y)\|_2^2\) where \(\gamma = \lambda \cdot T\).

  • The squared L2-norm gradient evaluated at x is given by \(2 \lambda \cdot At(w \cdot (A(x)-y))\).

Examples

>>> from pyunlocbox import functions
>>> f = functions.norm_l2()
>>> x = [1, 2, 3, 4]
>>> f.eval(x)
30
>>> f.prox(x, 1)
array([0.33333333, 0.66666667, 1.        , 1.33333333])
>>> f.grad(x)
array([2, 4, 6, 8])
class pyunlocbox.functions.norm_nuclear(**kwargs)[source]

Bases: pyunlocbox.functions.norm

Nuclear-norm (eval, prox).

See generic attributes descriptions of the pyunlocbox.functions.norm base class. Note that the constructor takes keyword-only parameters.

Notes

  • The nuclear-norm of the matrix x is given by \(\lambda \| x \|_* = \lambda \operatorname{trace} (\sqrt{x^* x}) = \lambda \sum_{i=1}^N |e_i|\) where e_i are the eigenvalues of x.

  • The nuclear-norm proximal operator evaluated at x is given by \(\operatorname{arg\,min}\limits_z \frac{1}{2} \|x-z\|_2^2 + \gamma \| x \|_*\) where \(\gamma = \lambda \cdot T\), which is a soft-thresholding of the eigenvalues.

Examples

>>> from pyunlocbox import functions
>>> f = functions.norm_nuclear()
>>> f.eval([[1, 2],[2, 3]])  
4.47213595...
>>> f.prox([[1, 2],[2, 3]], 1)
array([[0.89442719, 1.4472136 ],
       [1.4472136 , 2.34164079]])
class pyunlocbox.functions.norm_tv(dim=2, verbosity='LOW', **kwargs)[source]

Bases: pyunlocbox.functions.norm

TV-norm (eval, prox).

See generic attributes descriptions of the pyunlocbox.functions.norm base class. Note that the constructor takes keyword-only parameters.

Notes

TODO

See [BT09b] for details about the algorithm.

Examples

>>> import numpy as np
>>> from pyunlocbox import functions
>>> f = functions.norm_tv()
>>> x = np.arange(0, 16)
>>> x = x.reshape(4, 4)
>>> f.eval(x)  
    norm_tv evaluation: 5.210795e+01
52.10795063...
class pyunlocbox.functions.proj(**kwargs)[source]

Bases: pyunlocbox.functions.func

Base class which defines the attributes of the proj objects.

See generic attributes descriptions of the pyunlocbox.functions.func base class.

Notes

  • All indicator functions (projections) evaluate to zero by definition.

class pyunlocbox.functions.proj_positive(**kwargs)[source]

Bases: pyunlocbox.functions.proj

Projection on the positive octant (eval, prox).

This function is the indicator function \(i_S(z)\) of the set \(S = \left\{z \in \mathbb{R}^N \mid z \leq 0 \right\}\) that is zero if \(z\) is in the set and infinite otherwise.

See generic attributes descriptions of the pyunlocbox.functions.proj base class. Note that the constructor takes keyword-only parameters.

Notes

  • The evaluation of this function is zero.

Examples

>>> from pyunlocbox import functions
>>> f = functions.proj_positive()
>>> x = [-2.5, 1.5]
>>> f.eval(x)
0
>>> f.prox(x, 0)
array([0. , 1.5])
class pyunlocbox.functions.proj_spsd(**kwargs)[source]

Bases: pyunlocbox.functions.proj

Projection on symmetric positive semi-definite matrices (eval, prox).

This function is the indicator function \(i_S(M)\) of the set \(S = \left\{M \in \mathbb{R}^{N \times N} \mid M \succeq 0, M=M^T \right\}\) that is zero if \(M\) is in the set and infinite otherwise.

See generic attributes descriptions of the pyunlocbox.functions.proj base class. Note that the constructor takes keyword-only parameters.

Notes

  • The evaluation of this function is zero.

Examples

>>> from pyunlocbox import functions
>>> f = functions.proj_spsd()
>>> A = np.array([[0, -1] , [-1, 1]])
>>> A = (A + A.T) / 2  # Symmetrize the matrix.
>>> np.linalg.eig(A)[0]
array([-0.61803399,  1.61803399])
>>> f.eval(A)
0
>>> Aproj = f.prox(A, 0)
>>> np.linalg.eig(Aproj)[0]
array([0.        , 1.61803399])
class pyunlocbox.functions.proj_b2(epsilon=1, method='FISTA', **kwargs)[source]

Bases: pyunlocbox.functions.proj

Projection on the L2-ball (eval, prox).

This function is the indicator function \(i_S(z)\) of the set \(S= \left\{z \in \mathbb{R}^N \mid \|Az-y\|_2 \leq \epsilon \right\}\) that is zero if \(z\) is in the set and infinite otherwise.

See generic attributes descriptions of the pyunlocbox.functions.proj base class. Note that the constructor takes keyword-only parameters.

Parameters
epsilonfloat, optional

The radius of the ball. Default is 1.

method{‘FISTA’, ‘ISTA’}, optional

The method used to solve the problem. It can be ‘FISTA’ or ‘ISTA’. Default is ‘FISTA’.

See also

proj_lineq

use instead of epsilon=0

Notes

  • The tol parameter is defined as the tolerance for the projection on the L2-ball. The algorithm stops if \(\frac{\epsilon}{1-tol} \leq \|y-A(z)\|_2 \leq \frac{\epsilon}{1+tol}\).

  • The evaluation of this function is zero.

  • The L2-ball proximal operator evaluated at x is given by \(\operatorname{arg\,min}\limits_z \frac{1}{2} \|x-z\|_2^2 + i_S(z)\) which has an identical solution as \(\operatorname{arg\,min}\limits_z \|x-z\|_2^2\) such that \(\|A(z)-y\|_2 \leq \epsilon\). It is thus a projection of the vector x onto an L2-ball of diameter epsilon.

Examples

>>> from pyunlocbox import functions
>>> f = functions.proj_b2(y=[1, 1])
>>> x = [3, 3]
>>> f.eval(x)
0
>>> f.prox(x, 0)
array([1.70710678, 1.70710678])
class pyunlocbox.functions.proj_lineq(A=None, pinvA=None, **kwargs)[source]

Bases: pyunlocbox.functions.proj

Projection on the plane satisfying the linear equality Az = y (eval, prox).

This function is the indicator function \(i_S(z)\) of the set \(S = \left\{z \in \mathbb{R}^N \mid Az = y \right\}\) that is zero if \(z\) is in the set and infinite otherwise.

The proximal operator is \(\operatorname{arg\,min}_z \| z - x \|_2 \text{ s.t. } Az = y\).

See generic attributes descriptions of the pyunlocbox.functions.proj base class. Note that the constructor takes keyword-only parameters.

See also

proj_b2

quadratic case

Notes

  • A parameter pinvA, the pseudo-inverse of A, must be provided if the parameter A is provided as an operator/callable (not a matrix).

  • The evaluation of this function is zero.

Examples

>>> from pyunlocbox import functions
>>> import numpy as np
>>> x = np.array([0, 0])
>>> A = np.array([[1, 1]])
>>> pinvA = np.linalg.pinv(A)
>>> y = np.array([1])
>>> f = functions.proj_lineq(A=A, pinvA=pinvA, y=y)
>>> sol = f.prox(x, 0)
>>> sol
array([0.5, 0.5])
>>> np.abs(A.dot(sol) - y) < 1e-15
array([ True])
class pyunlocbox.functions.structured_sparsity(lambda_=1, groups=[[]], weights=[0], **kwargs)[source]

Bases: pyunlocbox.functions.func

Structured sparsity (eval, prox).

The structured sparsity term that is defined in the work of Jenatton et al. 2011 Proximal methods for hierarchical sparse coding.

\[\Omega(x) = \lambda \cdot \sum_{g \in G} w_g \cdot \|x_g\|_2\]

See generic attributes descriptions of the pyunlocbox.functions.func base class.

Parameters
lambda_float, optional

The scaling factor of the function that corresponds to \(\lambda\). Must be a non-negative number.

groups: list of lists of integers

Each element encodes the indices of the vector belonging to a single group. Corresponds to \(G\).

weightsarray_like

Weight associated to each group. Corresponds to \(w_g\). Must have the same length as \(G\).

Examples

>>> from pyunlocbox import functions
>>> groups = [[0, 1], [3, 2, 4]]
>>> weights = [2, 1]
>>> f = functions.structured_sparsity(10, groups, weights)
>>> x = [2, 2.5, -0.5, 0.3, 0.01]
>>> f.eval(x)
69.86305169905782
>>> f.prox(x, 0.1)
array([0.7506099 , 0.93826238, 0.        , 0.        , 0.        ])