Acceleration class hierarchy

Acceleration scheme object interface

class pyunlocbox.acceleration.accel[source]

Bases: object

Defines the acceleration scheme object interface.

This class defines the interface of an acceleration scheme object intended to be passed to a solver inheriting from pyunlocbox.solvers.solver. It is intended to be a base class for standard acceleration schemes which will implement the required methods. It can also be instantiated by user code and dynamically modified for rapid testing. This class also defines the generic attributes of all acceleration scheme objects.

post()[source]

Post-processing specific to the acceleration scheme.

Mainly used to delete references added during initialization so that the garbage collector can free the memory. Gets called when pyunlocbox.solvers.solve() finishes running.

pre(functions, x0)[source]

Pre-processing specific to the acceleration scheme.

Gets called when pyunlocbox.solvers.solve() starts running.

update_sol(solver, objective, niter)[source]

Update the solution point for the next iteration.

Parameters:

solver : pyunlocbox.solvers.solver

Solver on which to act.

objective : list of floats

List of evaluations of the objective function since the beginning of the iterative process.

niter : int

Current iteration number.

Returns:

array_like

Updated solution point.

update_step(solver, objective, niter)[source]

Update the step size for the next iteration.

Parameters:

solver : pyunlocbox.solvers.solver

Solver on which to act.

objective : list of floats

List of evaluations of the objective function since the beginning of the iterative process.

niter : int

Current iteration number.

Returns:

float

Updated step size.

Dummy acceleration scheme

class pyunlocbox.acceleration.dummy[source]

Bases: pyunlocbox.acceleration.accel

Dummy acceleration scheme.

Used by default in most of the solvers. It simply returns unaltered the step size and solution point it receives.

Backtracking from quadratic approximation

class pyunlocbox.acceleration.backtracking(eta=0.5, **kwargs)[source]

Bases: pyunlocbox.acceleration.dummy

Backtracking based on a local quadratic approximation of the the smooth part of the objective.

Parameters:

eta : float

A number between 0 and 1 representing the ratio of the geometric sequence formed by successive step sizes. In other words, it establishes the relation step_new = eta * step_old. Default is 0.5.

Notes

This is the backtracking strategy used in the original FISTA paper, [BT09a].

Examples

>>> from pyunlocbox import functions, solvers, acceleration
>>> import numpy as np
>>> y = [4, 5, 6, 7]
>>> x0 = np.zeros(len(y))
>>> f1 = functions.norm_l1(y=y, lambda_=1.0)
>>> f2 = functions.norm_l2(y=y, lambda_=0.8)
>>> accel = acceleration.backtracking()
>>> solver = solvers.forward_backward(accel=accel, step=10)
>>> ret = solvers.solve([f1, f2], x0, solver, atol=1e-32, rtol=None)
Solution found after 4 iterations:
    objective function f(sol) = 0.000000e+00
    stopping criterion: ATOL
>>> ret['sol']
array([ 4.,  5.,  6.,  7.])

FISTA acceleration scheme

class pyunlocbox.acceleration.fista(**kwargs)[source]

Bases: pyunlocbox.acceleration.dummy

Acceleration scheme for forward-backward solvers.

Notes

This is the acceleration scheme proposed in the original FISTA paper, [BT09a].

Examples

>>> from pyunlocbox import functions, solvers, acceleration
>>> import numpy as np
>>> y = [4, 5, 6, 7]
>>> x0 = np.zeros(len(y))
>>> f1 = functions.norm_l2(y=y)
>>> f2 = functions.dummy()
>>> accel=acceleration.fista()
>>> solver = solvers.forward_backward(accel=accel, step=0.5)
>>> ret = solvers.solve([f1, f2], x0, solver, atol=1e-5)
Solution found after 15 iterations:
    objective function f(sol) = 4.957288e-07
    stopping criterion: ATOL
>>> ret['sol']
array([ 4.0002509 ,  5.00031362,  6.00037635,  7.00043907])

FISTA acceleration scheme with backtracking

class pyunlocbox.acceleration.fista_backtracking(eta=0.5, **kwargs)[source]

Bases: pyunlocbox.acceleration.backtracking, pyunlocbox.acceleration.fista

Acceleration scheme with backtracking for forward-backward solvers.

Notes

This is the acceleration scheme and backtracking strategy proposed in the original FISTA paper, [BT09a].

Examples

>>> from pyunlocbox import functions, solvers, acceleration
>>> import numpy as np
>>> y = [4, 5, 6, 7]
>>> x0 = np.zeros(len(y))
>>> f1 = functions.norm_l2(y=y)
>>> f2 = functions.dummy()
>>> accel=acceleration.fista_backtracking()
>>> solver = solvers.forward_backward(accel=accel, step=0.5)
>>> ret = solvers.solve([f1, f2], x0, solver, atol=1e-5)
Solution found after 15 iterations:
    objective function f(sol) = 4.957288e-07
    stopping criterion: ATOL
>>> ret['sol']
array([ 4.0002509 ,  5.00031362,  6.00037635,  7.00043907])

Regularized Nonlinear Acceleration (RNA)

class pyunlocbox.acceleration.regularized_nonlinear(k=10, lambda_=1e-06, adaptive=True, dolinesearch=True, forcedecrease=True, **kwargs)[source]

Bases: pyunlocbox.acceleration.dummy

Regularized nonlinear acceleration (RNA) for gradient descent.

Parameters:

k : int, optional

Number of points to keep in the buffer for computing the extrapolation. (Default is 10.)

lambda_ : float or list of floats, optional

Regularization parameter in the acceleration scheme. The user can pass a list of candidates, and the acceleration algorithm will pick the one that provides the best extrapolation. (Default is 1e-6.)

adaptive : boolean, optional

If adaptive = True and the user has not provided a list of regularization parameters, the acceleration algorithm will assemble a grid of possible regularization parameters based on the SVD of the Gram matrix of vectors of differences in the extrapolation buffer. If adaptive = False, the algorithm will simply try to use the value(s) given in lambda_. (Default is True.)

dolinesearch : boolean, optional

If dolinesearch = True, the acceleration scheme will try to return a point in the line segment between the current extrapolation and the previous one that provides a decrease in the value of the objective function. If dolinesearch = False, the algorithm simply returns the current extrapolation. (Default is True.)

forcedecrease : boolean, optional

If forcedecrese = True and we obtain a bad extrapolation, the algorithm returns the unchanged solution produced by the solver. If forcedecrease = False, the algorithm returns the new extrapolation no matter what. (Default is True.)

Notes

This is the acceleration scheme proposed in [SdB16].

See also Damien Scieur’s repository for the Matlab version that inspired this implementation.

Examples

>>> from pyunlocbox import functions, solvers, acceleration
>>> import numpy as np
>>> dim = 25;
>>> np.random.seed(0)
>>> xstar = np.random.rand(dim) # True solution
>>> x0 = np.random.rand(dim)
>>> x0 = xstar + 5.*(x0 - xstar) / np.linalg.norm(x0 - xstar)
>>> A = np.random.rand(dim, dim)
>>> step = 1/np.linalg.norm(np.dot(A.T, A))
>>> f = functions.norm_l2(lambda_=0.5, A=A, y=np.dot(A, xstar))
>>> fd = functions.dummy()
>>> accel = acceleration.regularized_nonlinear(k=5)
>>> solver = solvers.gradient_descent(step=step, accel=accel)
>>> params = {'rtol':0, 'maxit':200, 'verbosity':'NONE'}
>>> ret = solvers.solve([f, fd], x0, solver, **params)
>>> pctdiff = 100*np.sum((xstar - ret['sol'])**2)/np.sum(xstar**2)
>>> print('Difference: {0:.1f}%'.format(pctdiff))
Difference: 1.3%
lambda_