regpy.functionals.base

Classes

AbstractFunctional

An abstract functional that can be called on a vector space to get the corresponding

Functional

Base class for implementation of functionals. Subclasses should at least implement the

LinearFunctional

Linear functionals

SquaredNorm

Functionals of the form

LinearCombination

Linear combination of functionals.

VerticalShift

Shifting a functional by some offset. Should not be used directly but rather by adding some scalar to the functional.

HorizontalShiftDilation

Implements a horizontal shift and/or a horizontal translation of the graph of a functional \(F\), i.e. replaces

Composed

Composition of an operator with a functional \(F\circ O\). This should not be called

FunctionalOnDirectSum

Helper to define Functionals with respective prox-operators on sum spaces (vecsps.DirectSum objects).

Functions

as_functional(func, vecsp)

Convert func to Functional instance on vecsp.

Module Contents

class regpy.functionals.base.AbstractFunctional(name)[source]

Bases: AbstractFunctionalBase

An abstract functional that can be called on a vector space to get the corresponding concrete implementation.

AbstractFunctionals provides two kinds of functionality:

  • A decorator method register(vecsp_type) that can be used to declare some class or function as the concrete implementation of this abstract functional for vector spaces of type vecsp_type or subclasses thereof, e.g.:

  • AbstractFunctionals are callable. Calling them on a vector space and arbitrary optional keyword arguments finds the corresponding concrete regpy.functionals.Functional among all registered implementations. If there are implementations for multiple base classes of the vector space type, the most specific one will be chosen. The chosen implementation will then be called with the vector space and the keyword arguments, and the result will be returned.

@TV.register(vecsps.UniformGridFcts)
class TVUniformGridFcts(HilbertSpace):
    ...

If called without a vector space as positional argument, it returns a new abstract functional with all passed keyword arguments remembered as defaults.

Parameters:

name (str) – A name for this abstract functional. Currently, this is only used in error messages, when no implementation was found for some vector space.

name
args
register(vecsp_type, impl=None)[source]

Either registers a new implementation on a specific regpy.vecsps.VectorSpaceBase for a given Abstract functional or returns as decorator that can output any implementation option for a given vector space.

Parameters:
  • vecsp_type (regpy.vecsps.VectorSpaceBase) – Vector Space on which the functional should be registered.

  • impl (regpy.functionals.Functional, optional) – The explicit implementation to be used for that Vector Space, by default None

Returns:

None or decorator – Either nothing or map that can output any of the registered implementations for a specific vector space.

Return type:

None or map

class regpy.functionals.base.Functional(domain, h_domain=None, linear=False, convexity_param=0.0, Lipschitz=inf, separable=False, convex=True, dom_l=None, dom_u=None, conj_dom_l=None, conj_dom_u=None, methods=set(), conj_methods=set(), is_data_func=False)[source]
Base class for implementation of functionals. Subclasses should at least implement the

_eval : evaluating the functional

and

_subgradient or _linearize : returning a subgradient at x.

The evaluation of a specific functional on some element of the domain can be done by simply calling the functional on that element.

Functionals can be added by taking LinearCombination of them. The domain has to be the same for each functional.

They can also be multiplied by scalars or vector of their respective domain`or multiplied by `regpy.operators.Operator. This leads to a functional that is composed with the operator \(F\circ O\) where \(F\) is the functional and \(O\) some operator. Multiplying by a scalar results in a composition with the PtwMultiplication operator.

Parameters:
  • domain (regpy.vecsps.VectorSpaceBase) – The underlying vector space for the function space on which it is defined.

  • h_domain (regpy.hilbert.HilbertSpace (default: None)) – The underlying Hilbert space. The proximal mapping, the parameter of strong convexity, and the Lipschitz constant are defined with respect to this Hilbert space. In the default case L2(domain) is used.

  • convex (bool [default: True]) – If true, the functional should be convex.

  • linear (bool [default: False]) – If true, the functional should be linear.

  • separable (bool [default: False]) – If true, the functional should be the sum of functionals acting on only one component of the input vector. In this case, the parameters

  • dom_u (self.domain [default:None]) – should not be None, and they should specify the essential domain of the functional by \(\{x in domain: dom_l<=x<=dom_u}\), and the essential domain of the conjugate functional (which is then also separable) by \(\{xstar in domain: conj_dom_l<=xstar <= conj_dom_u\}\). In case of open domains, the boundaries should be shifted in the order of machine precision.

  • dom_l (self.domain [default:None]) – should not be None, and they should specify the essential domain of the functional by \(\{x in domain: dom_l<=x<=dom_u}\), and the essential domain of the conjugate functional (which is then also separable) by \(\{xstar in domain: conj_dom_l<=xstar <= conj_dom_u\}\). In case of open domains, the boundaries should be shifted in the order of machine precision.

  • conj_dom_u (self.domain [default:None]) – should not be None, and they should specify the essential domain of the functional by \(\{x in domain: dom_l<=x<=dom_u}\), and the essential domain of the conjugate functional (which is then also separable) by \(\{xstar in domain: conj_dom_l<=xstar <= conj_dom_u\}\). In case of open domains, the boundaries should be shifted in the order of machine precision.

  • conj_dom_l (self.domain [default:None]) – should not be None, and they should specify the essential domain of the functional by \(\{x in domain: dom_l<=x<=dom_u}\), and the essential domain of the conjugate functional (which is then also separable) by \(\{xstar in domain: conj_dom_l<=xstar <= conj_dom_u\}\). In case of open domains, the boundaries should be shifted in the order of machine precision.

  • convexity_param (float [default: 0]) – parameter of strong convexity of the functional. 0 if the functional is not strongly convex.

  • Lipschitz (float [default: math.inf]) – Lipschitz continuity constant of the gradient. math.inf the gradient is not Lipschitz continuous.

  • methods (set [default: set()]) – names of the methods implemented by a given Functional instance. Subset of {‘eval’, ‘subgradient’, ‘hessian’, ‘proximal’, ‘dist_subdiff’}

  • [default (conj_methods set) – names of the methods implemented by the conjugate of a given Functional instance. Subset of {‘eval’, ‘subgradient’, ‘hessian’, ‘proximal’, ‘dist_subdiff’}

log
domain

The underlying vector space.

h_domain

The underlying Hilbert space.

convexity_param

parameter of strong convexity of the functional.

Lipschitz

Lipschitz continuity constant of the gradient.

linear = False

boolean indicating if the functional is linear

separable = False

boolean indicating if the functional is separable.

convex = True

boolean indicating if the functional is convex.

is_data_func = False
linearize(x)[source]

Bounds the functional from below by a linear functional at x given by the value at that point and a subgradient v such that

\[F(x+ h) \geq F(x) + vdot(v,h) for all h\]

Requires the implementation of either _subgradient or _linearize.

Parameters:

x (in self.domain) – Element at which will be linearized

Returns:

  • y – Value of \(F(x)\).

  • grad (in self.domain) – Subgradient of \(F\) at \(x\).

subgradient(x)[source]

Returns a subgradient \(\xi\) of the functional at x characterized by

\[F(y) \geq F(x) + vdot(\xi,y-x) for all y\]

Requires the implementation of either _subgradient or _linearize.

Parameters:

x (in self.domain) – Element at which will be linearized

Returns:

grad – subgradient of \(F\) at \(x\).

Return type:

in self.domain

dist_subdiff(vstar, x)[source]

Returns the distance of a vector \(v^*\) to the subdifferential \(\partial F(x)\) at x with respect to the dual norm. Needs to be re-implemented for functionals which are not Gateaux differentiable.

Parameters:
  • vstar (in self.domain) – Vector of which the distance is to be determined.

  • x (in self.domain) – Point at which the subdifferential is evaluated

hessian(x, recursion_safeguard=False)[source]

The hessian of the functional at x as an regpy.operators.Operator mapping form the functionals domain to it self. It is defined by

\[F(x+h) = F(x) + (\nabla F)(x)^T h + \frac{1}{2} h^T Hess F(x) h + \mathcal{o}(\|h\|^2)\]

Require either the implementation of _hessian or of _hessian_conj and _subgradient

Parameters:

x (self.domain) – Point in domain at which to compute the hessian.

Returns:

`h` – Hessian operator at the point x.

Return type:

regpy.operators.Operator

conj_subgradient(xstar)[source]

Gradient of the conjugate functional. Should not be called directly, but via self.conj.subgradient. Requires the implementation of _conj_subgradient.

conj_hessian(xstar, recursion_safeguard=False)[source]

The hessian of the functional. Should not be called directly, but via self.conj.hessian.

conj_linearize(xstar)[source]

Linearizes the conjugate functional \(F^*\). Should not be called directly, but via self.conj.linearize

proximal(x, tau, recursion_safeguard=False, **proximal_par)[source]

Proximal operator

\[\mathrm{prox}_{\tau F}(x)=\arg \min _{v\in {\mathcal {X}}}(F(v)+{\frac{1}{2\tau}}\Vert v-x\Vert_{\mathcal {X}}^{2}).\]

Requires either an implementation of _proximal or of _subgradient and _conj_proximal.

Parameters:
  • x (array-like) – Vector in the respective domain. Point at which to compute proximal.

  • tau (scalar) – Regularization parameter for the proximal.

Returns:

proximal – the computed proximal at \(x\) with parameter \(\tau\).

Return type:

self.domain

conj_proximal(xstar, tau, recursion_safeguard=False, **proximal_par)[source]

Proximal operator of conjugate functional. Should not be called directly, but via self.conj.proximal

as_data_func(data)[source]
shift(v=None, data_shift=None)[source]

Returns the functional \(x\mapsto F(x-v-data)\)

dilation(a)[source]

Returns the functional \(x\mapsto F(ax)\)

property methods

Set of strings of the names of the methods are available for a given Functional instance.

conj()[source]

For linear operators, this is the adjoint as a linear regpy.operators.Operator instance. Will only be computed on demand and saved for subsequent invocations.

Returns:

The adjoint as an regpy.operators.Operator instance.

Return type:

Adjoint

class regpy.functionals.base.LinearFunctional(gradient, domain=None, h_domain=None, gradient_in_dual_space=False)[source]

Bases: Functional

Linear functionals Linear functional given by

\[F(x) = \langel a, x\rangle\]

The operators __add__ , __iadd__ , __mul__ , __imul__ with LinearFunctionals and scalars as other arguments, rsp, are overwritten to yield the expected LinearFunctionals.

Parameters:
  • gradient (domain) – The gradient of the linear functional. \(a=gradient\) if gradient_in_dual_space == True

  • domain (regpy.vecsps.VectorSpaceBase, optional) – The regpy.vecsps.VectorSpaceBase on which the functional is defined

  • h_domain (regpy.hilbert.HilbertSpace (default: L2(domain))) – Hilbert space for proximity operator

  • gradient_in_dual_space (bool (default: False)) – If false, the argument gradient is considered as an element of the primal space, and \(a = h_domain.gram(gradient)\).

property gradient
dilation(a)[source]

Returns the functional \(x\mapsto F(ax)\)

shift(v=None, data_shift=None)[source]

Returns the functional \(x\mapsto F(x-v-data)\)

class regpy.functionals.base.SquaredNorm(h_space, a=1.0, b=None, c=0.0, shift=None, data=None)[source]

Bases: Functional

Functionals of the form

\[\mathcal{F}(x) = \frac{a}{2}\|x\|_X^2 +\langle b,x\rangle_X + c\]

Here the linear term represents an inner product in the Hilbert space, not a pairing with the dual space.

The operators __add__ , __iadd__ , __mul__ , __imul__ with SquaredNorms, LinearFunctionals and scalars as other arguments are overwritten to yield the expected SquaredNorms.

Parameters:
  • domain (regpy.vecsps.VectorSpaceBase) – The uncerlying vector space for the function space on which it is defined.

  • h_sapce (regpy.hilbert.HilbertSpace (default: None)) – The underlying Hilbert space.

  • a (float [default:1]) – coefficient of quadratic term

  • b (h_space.domain [default:None]) – coefficient of linear term. In the default case it is 0.

  • c (float [default: 0]) – constant term

  • shift (h_space.domain [default:None]) – If not None, then we must have b is None and c==0. In this case the functional is initialized as \(\mathcal{F}(x) = \frac{a}{2}\|x-shift-data\|^2\).

  • data (h_space.domain [default:None]) – If not None, then we must have b is None and c==0. In this case the functional is initialized as \(\mathcal{F}(x) = \frac{a}{2}\|x-shift-data\|^2\).

gram
property data
a()[source]
b()[source]
c()[source]
as_data_func(data)[source]
dilation(dil)[source]

Returns the functional \(x\mapsto F(ax)\)

shift(v=None, data_shift=None)[source]

Returns the functional \(x\mapsto F(x-v-data)\)

class regpy.functionals.base.LinearCombination(*args)[source]

Bases: Functional

Linear combination of functionals.

Parameters:

*args ((scalar, regpy.functionals.Functional) or regpy.functionals.Functional) – List of coefficients and functionals to be taken as linear combinations.

coeffs = []

List of all coefficients

funcs = []

List of all functionals.

linear_table = []
dist_subdiff(vstar, x, **kwargs)[source]

Returns the distance of a vector \(v^*\) to the subdifferential \(\partial F(x)\) at x with respect to the dual norm. Needs to be re-implemented for functionals which are not Gateaux differentiable.

Parameters:
  • vstar (in self.domain) – Vector of which the distance is to be determined.

  • x (in self.domain) – Point at which the subdifferential is evaluated

class regpy.functionals.base.VerticalShift(func, offset)[source]

Bases: Functional

Shifting a functional by some offset. Should not be used directly but rather by adding some scalar to the functional.

Parameters:
  • func (regpy.functionals.Functional) – Functional to be offset.

  • offset (scalar) – Reals offset added to the evaluation of the functional.

func

Functional to be offset.

offset

Offset added to the evaluation of the functional.

dist_subdiff(vstar, x, **kwargs)[source]

Returns the distance of a vector \(v^*\) to the subdifferential \(\partial F(x)\) at x with respect to the dual norm. Needs to be re-implemented for functionals which are not Gateaux differentiable.

Parameters:
  • vstar (in self.domain) – Vector of which the distance is to be determined.

  • x (in self.domain) – Point at which the subdifferential is evaluated

class regpy.functionals.base.HorizontalShiftDilation(func, dilation=1.0, shift=None, data=None)[source]

Bases: Functional

Implements a horizontal shift and/or a horizontal translation of the graph of a functional \(F\), i.e. replaces \(F(x)\) by \(F(dilation(x-shift))\)

Parameters:
  • func (Functional) – The functional to be shifted and dilated.

  • dilation (float [default: 1]) – Dilation factor.

  • shift (self.domain or scalar or None [default: None]) – Shift vector. The default case (None) yields the same results as shift=0, but no zero-additions are performed.

func
dilation = 1.0

Returns the functional \(x\mapsto F(ax)\)

shift_val()[source]
recompute_cutoff()[source]
property data
dist_subdiff(vstar, x, **kwargs)[source]

Returns the distance of a vector \(v^*\) to the subdifferential \(\partial F(x)\) at x with respect to the dual norm. Needs to be re-implemented for functionals which are not Gateaux differentiable.

Parameters:
  • vstar (in self.domain) – Vector of which the distance is to be determined.

  • x (in self.domain) – Point at which the subdifferential is evaluated

class regpy.functionals.base.Composed(func, op, op_norm=inf, op_lower_bound=0, compute_op_norm=False, norm_kwargs={}, methods=None, conj_methods=None)[source]

Bases: Functional

Composition of an operator with a functional \(F\circ O\). This should not be called directly but rather used by multiplying the Functional object with an regpy.operators.Operator.

Parameters:
  • func (regpy.functionals.Functional) – Functional to be composed with.

  • op (regpy.operators.Operator) – Operator to be composed with.

  • op_norm (float [default: inf]) – Norm of the operator. Used only to define self.Lipschitz

  • op_lower_bound (float) – Lower bound of operator: |op(f)|geq op_lower_bound * |f| Used only to define self.convexity_param

  • compute_op_norm (boolean) – If true the op norm will be computed using the norm method of the operator. Will only be computed if op_norm is default value inf.

  • norm_kwargs (dict) – possible arguments passed to the operator norm computation.

func

Functional that is composed with an Operator.

op

Operator composed that is composed with a functional.

class regpy.functionals.base.FunctionalOnDirectSum(funcs, domain=None)[source]

Bases: Functional

Helper to define Functionals with respective prox-operators on sum spaces (vecsps.DirectSum objects). The functionals are given as a list of the functionals on the summands of the sum space.

\[F(x_1,... x_n) = \sum_{j=1}^n F_j(x_j)\]
Parameters:
length

Number of the summands in the direct sum domain.

funcs

List of the functionals on each summand of the direct sum domain.

dist_subdiff(vstar, x)[source]

Returns the distance of a vector \(v^*\) to the subdifferential \(\partial F(x)\) at x with respect to the dual norm. Needs to be re-implemented for functionals which are not Gateaux differentiable.

Parameters:
  • vstar (in self.domain) – Vector of which the distance is to be determined.

  • x (in self.domain) – Point at which the subdifferential is evaluated

regpy.functionals.base.as_functional(func, vecsp)[source]

Convert func to Functional instance on vecsp.

  • If func is a HilbertSpace then it generated the SquaredNorm.

  • If func is an Operator, it’s wrapped in a GramHilbertSpace and then SquaredNorm functional.

  • If func is callable, e.g. an hilbert.AbstractSpace or AbstractFunctional, it is called on vecsp to construct the concrete functional or Hilbert space. In the later case the functional will be the SquaredNorm

Parameters:
  • func (Functional or HilbertSapce or regpy.operators.Operator or callable) – Functional or object from which to construct the Functional.

  • vecsp (regpy.vecsps.VectorSpaceBase) – Underlying vector space for the functional.

Returns:

Constructed Functional on the underlying vectorspace.

Return type:

Functional