/
functions.py
896 lines (735 loc) · 27.8 KB
/
functions.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
# -*- coding: utf-8 -*-
r"""
This module implements function objects which are then passed to solvers. The
:class:`func` base class defines the interface whereas specialised classes who
inherit from it implement the methods. These classes include :
* :class:`dummy`: A dummy function object which returns 0 for the
:meth:`_eval`, :meth:`_prox` and :meth:`_grad` methods.
* :class:`norm`: Norm operators base class.
* :class:`norm_l1`: L1-norm who implements the :meth:`_eval` and
:meth:`_prox` methods.
* :class:`norm_l2`: L2-norm who implements the :meth:`_eval`, :meth:`_prox`
and :meth:`_grad` methods.
* :class:`norm_nuclear`: nuclear-norm who implements the :meth:`_eval` and
:meth:`_prox` methods.
* :class:`norm_tv`: TV-norm who implements the :meth:`_eval` and
:meth:`_prox` methods.
* :class:`proj`: Projection operators base class.
* :class:`proj_b2`: Projection on the L2-ball who implements the
:meth:`_eval` and :meth:`_prox` methods.
"""
from time import time
import numpy as np
from copy import deepcopy
from scipy.optimize import minimize
from pyunlocbox import operators as op
def _soft_threshold(z, T, handle_complex=True):
r"""
Return the soft thresholded signal.
Parameters
----------
z : array_like
Input signal (real or complex).
T : float or array_like
Threshold on the absolute value of `z`. There could be either a single
threshold for the entire signal `z` or one threshold per dimension.
Useful when you use weighted norms.
handle_complex : bool
Indicate that we should handle the thresholding of complex numbers,
which may be slower. Default is True.
Returns
-------
sz : ndarray
Soft thresholded signal.
Examples
--------
>>> import pyunlocbox
>>> pyunlocbox.functions._soft_threshold([-2, -1, 0, 1, 2], 1)
array([-1, 0, 0, 0, 1])
"""
sz = np.maximum(np.abs(z) - T, 0)
if not handle_complex:
# This soft thresholding method only supports real signal.
sz[:] = np.sign(z) * sz
else:
# This soft thresholding method supports complex complex signal.
# Transform to float to avoid integer division.
# In our case 0 divided by 0 should be 0, not NaN, and is not an error.
# It corresponds to 0 thresholded by 0, which is 0.
old_err_state = np.seterr(invalid='ignore')
sz[:] = np.nan_to_num(1. * sz / (sz + T) * z)
np.seterr(**old_err_state)
return sz
class func(object):
r"""
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 :func:`pyunlocbox.solvers.solve` solving
function.
Parameters
----------
y : array_like, optional
Measurements. Default is 0.
A : function or ndarray, optional
The forward operator. Default is the identity, :math:`A(x)=x`. If `A`
is an ``ndarray``, it will be converted to the operator form.
At : function 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`,
:math:`At(x)=A(x)`.
tight : bool, optional
``True`` if `A` is a tight frame, ``False`` otherwise. Default is
``True``.
nu : float, optional
Bound on the norm of the operator `A`, i.e. :math:`\|A(x)\|^2 \leq \nu
\|x\|^2`. Default is 1.
tol : float, 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.
maxit : int, 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 :
>>> import pyunlocbox
>>> f = pyunlocbox.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']
"""
def __init__(self, y=0, A=None, At=None, tight=True, nu=1, tol=1e-3,
maxit=200, **kwargs):
if callable(y):
self.y = lambda: np.asarray(y())
else:
self.y = lambda: np.asarray(y)
if A is None:
self.A = lambda x: x
else:
if callable(A):
self.A = A
else:
# Transform matrix form to operator form.
self.A = lambda x: np.dot(A, x)
if At is None:
if A is None:
self.At = lambda x: x
elif callable(A):
self.At = A
else:
self.At = lambda x: np.dot(np.transpose(A), x)
else:
if callable(At):
self.At = At
else:
self.At = lambda x: np.dot(At, x)
self.tight = tight
self.nu = nu
self.tol = tol
self.maxit = maxit
# Should be initialized if called alone, updated by solve().
self.verbosity = 'NONE'
def eval(self, x):
r"""
Function evaluation.
Parameters
----------
x : array_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
-------
z : float
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 :func:`pyunlocbox.solvers.solve` solving
function to evaluate the objective function. Each function class
should therefore define it.
"""
sol = self._eval(np.asarray(x))
if self.verbosity in ['LOW', 'HIGH']:
name = self.__class__.__name__
print(' {} evaluation: {:e}'.format(name, sol))
return sol
def _eval(self, x):
raise NotImplementedError("Class user should define this method.")
def prox(self, x, T):
r"""
Function proximal operator.
Parameters
----------
x : array_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.
T : float
The regularization parameter.
Returns
-------
z : ndarray
The proximal operator evaluated for each column of `x`.
Notes
-----
This method is required by some solvers.
The proximal operator is defined by
:math:`\operatorname{prox}_{\gamma f}(x) = \operatorname{arg\,min}
\limits_z \frac{1}{2} \|x-z\|_2^2 + \gamma f(z)`
"""
return self._prox(np.asarray(x), T)
def _prox(self, x, T):
raise NotImplementedError("Class user should define this method.")
def grad(self, x):
r"""
Function gradient.
Parameters
----------
x : array_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
-------
z : ndarray
The objective function gradient evaluated for each column of `x`.
Notes
-----
This method is required by some solvers.
"""
return self._grad(np.asarray(x))
def _grad(self, x):
raise NotImplementedError("Class user should define this method.")
def cap(self, x):
r"""
Test the capabilities of the function object.
Parameters
----------
x : array_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
-------
cap : list of string
A list of capabilities ('EVAL', 'GRAD', 'PROX').
"""
tmp = self.verbosity
self.verbosity = 'NONE'
cap = ['EVAL', 'GRAD', 'PROX']
try:
self.eval(x)
except NotImplementedError:
cap.remove('EVAL')
try:
self.grad(x)
except NotImplementedError:
cap.remove('GRAD')
try:
self.prox(x, 1)
except NotImplementedError:
cap.remove('PROX')
self.verbosity = tmp
return cap
class dummy(func):
r"""
Dummy function object.
This can be used as a second function object when there is only one
function to minimize. It always evaluates as 0.
Examples
--------
>>> import pyunlocbox
>>> f = pyunlocbox.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.])
"""
def __init__(self, **kwargs):
# Constructor takes keyword-only parameters to prevent user errors.
super(dummy, self).__init__(**kwargs)
def _eval(self, x):
return 0
def _prox(self, x, T):
return x
def _grad(self, x):
return np.zeros(np.shape(x))
class norm(func):
r"""
Base class which defines the attributes of the `norm` objects.
See generic attributes descriptions of the
:class:`pyunlocbox.functions.func` base class.
Parameters
----------
lambda_ : float, optional
Regularization parameter :math:`\lambda`. Default is 1.
w : array_like, optional
Weights for a weighted norm. Default is 1.
"""
def __init__(self, lambda_=1, w=1, **kwargs):
super(norm, self).__init__(**kwargs)
self.lambda_ = lambda_
self.w = np.asarray(w)
class norm_l1(norm):
r"""
L1-norm function object.
See generic attributes descriptions of the
:class:`pyunlocbox.functions.norm` base class. Note that the constructor
takes keyword-only parameters.
Notes
-----
* The L1-norm of the vector `x` is given by
:math:`\lambda \|w \cdot (A(x)-y)\|_1`.
* The L1-norm proximal operator evaluated at `x` is given by
:math:`\operatorname{arg\,min}\limits_z \frac{1}{2} \|x-z\|_2^2 + \gamma
\|w \cdot (A(z)-y)\|_1` where :math:`\gamma = \lambda \cdot T`. This is
simply a soft thresholding.
Examples
--------
>>> import pyunlocbox
>>> f = pyunlocbox.functions.norm_l1()
>>> f.eval([1, 2, 3, 4])
10
>>> f.prox([1, 2, 3, 4], 1)
array([0, 1, 2, 3])
"""
def __init__(self, **kwargs):
# Constructor takes keyword-only parameters to prevent user errors.
super(norm_l1, self).__init__(**kwargs)
def _eval(self, x):
sol = self.A(x) - self.y()
return self.lambda_ * np.sum(np.abs(self.w * sol))
def _prox(self, x, T):
# Gamma is T in the matlab UNLocBox implementation.
gamma = self.lambda_ * T
if self.tight:
# Nati: I've checked this code the use of 'y' seems correct
sol = self.A(x) - self.y()
sol[:] = _soft_threshold(sol, gamma * self.nu * self.w) - sol
sol[:] = x + self.At(sol) / self.nu
else:
raise NotImplementedError('Not implemented for non tight frame.')
return sol
class norm_l2(norm):
r"""
L2-norm function object.
See generic attributes descriptions of the
:class:`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
:math:`\lambda \|w \cdot (A(x)-y)\|_2^2`.
* The squared L2-norm proximal operator evaluated at `x` is given by
:math:`\operatorname{arg\,min}\limits_z \frac{1}{2} \|x-z\|_2^2 + \gamma
\|w \cdot (A(z)-y)\|_2^2` where :math:`\gamma = \lambda \cdot T`.
* The squared L2-norm gradient evaluated at `x` is given by
:math:`2 \lambda \cdot At(w \cdot (A(x)-y))`.
Examples
--------
>>> import pyunlocbox
>>> f = pyunlocbox.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])
"""
def __init__(self, **kwargs):
# Constructor takes keyword-only parameters to prevent user errors.
super(norm_l2, self).__init__(**kwargs)
def _eval(self, x):
sol = self.A(x) - self.y()
return self.lambda_ * np.sum((self.w * sol)**2)
def _prox(self, x, T):
# Gamma is T in the matlab UNLocBox implementation.
gamma = self.lambda_ * T
if self.tight:
sol = x + 2. * gamma * self.At(self.y() * self.w**2)
sol /= 1. + 2. * gamma * self.nu * self.w**2
else:
res = minimize(fun=lambda z: 0.5 * np.sum((z - x)**2) + gamma *
np.sum((self.w * (self.A(z) - self.y()))**2),
x0=x,
method='BFGS',
jac=lambda z: z - x + 2. * gamma *
self.At((self.w**2) * (self.A(z) - self.y())))
if res.success:
sol = res.x
else:
raise RuntimeError('norm_l2.prox: ' + res.message)
return sol
def _grad(self, x):
sol = self.A(x) - self.y()
return 2 * self.lambda_ * self.At((self.w**2) * sol)
class norm_nuclear(norm):
r"""
Nuclear-norm function object.
See generic attributes descriptions of the
:class:`pyunlocbox.functions.norm` base class. Note that the constructor
takes keyword-only parameters.
Notes
-----
* The nuclear-norm of the matrix `x` is given by
:math:`\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
:math:`\operatorname{arg\,min}\limits_z \frac{1}{2} \|x-z\|_2^2 + \gamma
\| x \|_*` where :math:`\gamma = \lambda \cdot T`, which is a
soft-thresholding of the eigenvalues.
Examples
--------
>>> import pyunlocbox
>>> f = pyunlocbox.functions.norm_nuclear()
>>> f.eval([[1, 2],[2, 3]]) # doctest:+ELLIPSIS
4.47213595...
>>> f.prox([[1, 2],[2, 3]], 1)
array([[ 0.89442719, 1.4472136 ],
[ 1.4472136 , 2.34164079]])
"""
def __init__(self, **kwargs):
# Constructor takes keyword-only parameters to prevent user errors.
super(norm_nuclear, self).__init__(**kwargs)
def _eval(self, x):
# TODO: take care of sparse matrices.
_, s, _ = np.linalg.svd(x)
return self.lambda_ * np.sum(np.abs(s))
def _prox(self, x, T):
# Gamma is T in the matlab UNLocBox implementation.
gamma = self.lambda_ * T
# TODO: take care of sparse matrices.
U, s, V = np.linalg.svd(x)
s = _soft_threshold(s, gamma)
S = np.diag(s)
return np.dot(U, np.dot(S, V))
class norm_tv(norm):
r"""
TV Norm function object.
See generic attributes descriptions of the
:class:`pyunlocbox.functions.norm` base class. Note that the constructor
takes keyword-only parameters.
Notes
-----
TODO
See :cite:`beck2009fastTV` for details about the algorithm.
Examples
--------
>>> import pyunlocbox
>>> import numpy as np
>>> f = pyunlocbox.functions.norm_tv()
>>> x = np.arange(0, 16)
>>> x = x.reshape(4, 4)
>>> f.eval(x) # doctest:+ELLIPSIS
norm_tv evaluation: 5.210795e+01
52.10795063...
"""
def __init__(self, dim=2, verbosity='LOW', **kwargs):
super(norm_tv, self).__init__(**kwargs)
self.kwargs = kwargs
self.dim = dim
self.verbosity = verbosity
def _eval(self, x):
if self.dim >= 2:
y = 0
grads = []
grads = op.grad(x, dim=self.dim, **self.kwargs)
for g in grads:
y += np.power(abs(g), 2)
y = np.sqrt(y)
return np.sum(y)
if self.dim == 1:
dx = op.grad(x, dim=self.dim, **self.kwargs)
y = np.sum(np.abs(dx), axis=0)
return np.sum(y)
def _prox(self, x, T):
# Time counter
t_init = time()
tol = self.tol
maxit = self.maxit
# TODO implement test_gamma
# Initialization
sol = x
if self.dim == 1:
r = op.grad(x * 0, dim=self.dim, **self.kwargs)
rr = deepcopy(r)
elif self.dim == 2:
r, s = op.grad(x * 0, dim=self.dim, **self.kwargs)
rr, ss = deepcopy(r), deepcopy(s)
elif self.dim == 3:
r, s, k = op.grad(x * 0, dim=self.dim, **self.kwargs)
rr, ss, kk = deepcopy(r), deepcopy(s), deepcopy(k)
elif self.dim == 4:
r, s, k, u = op.grad(x * 0, dim=self.dim, **self.kwargs)
rr, ss, kk, uu = deepcopy(r), deepcopy(s), deepcopy(k), deepcopy(u)
if self.dim >= 1:
pold = r
if self.dim >= 2:
qold = s
if self.dim >= 3:
kold = k
if self.dim >= 4:
uold = u
told, prev_obj = 1., 0.
# Initialization for weights
if self.dim >= 1:
try:
wx = self.kwargs["wx"]
except (KeyError, TypeError):
wx = 1.
if self.dim >= 2:
try:
wy = self.kwargs["wy"]
except (KeyError, TypeError):
wy = 1.
if self.dim >= 3:
try:
wz = self.kwargs["wz"]
except (KeyError, TypeError):
wz = 1.
if self.dim >= 4:
try:
wt = self.kwargs["wt"]
except (KeyError, TypeError):
wt = 1.
if self.dim == 1:
mt = wx
elif self.dim == 2:
mt = np.maximum(wx, wy)
elif self.dim == 3:
mt = np.maximum(wx, np.maximum(wy, wz))
elif self.dim == 4:
mt = np.maximum(np.maximum(wx, wy), np.maximum(wz, wt))
if self.verbosity in ['LOW', 'HIGH', 'ALL']:
print("Proximal TV Operator")
iter = 0
while iter <= maxit:
# Current Solution
if self.dim == 1:
sol = x - T * op.div(rr, **self.kwargs)
elif self.dim == 2:
sol = x - T * op.div(rr, ss, **self.kwargs)
elif self.dim == 3:
sol = x - T * op.div(rr, ss, kk, **self.kwargs)
elif self.dim == 4:
sol = x - T * op.div(rr, ss, kk, uu, **self.kwargs)
# Objective function value
obj = 0.5 * np.power(np.linalg.norm(x[:] - sol[:]), 2) + \
T * np.sum(self._eval(sol), axis=0)
rel_obj = np.abs(obj - prev_obj) / obj
prev_obj = obj
if self.verbosity in ['HIGH', 'ALL']:
print("Iter: ", iter, " obj = ", obj, " rel_obj = ", rel_obj)
# Stopping criterion
if rel_obj < tol:
crit = "TOL_EPS"
break
# Update divergence vectors and project
if self.dim == 1:
dx = op.grad(sol, dim=self.dim, **self.kwargs)
r -= 1. / (4 * T * mt**2) * dx
weights = np.maximum(1, np.abs(r))
elif self.dim == 2:
dx, dy = op.grad(sol, dim=self.dim, **self.kwargs)
r -= (1. / (8. * T * mt**2.)) * dx
s -= (1. / (8. * T * mt**2.)) * dy
weights = np.maximum(1, np.sqrt(np.power(np.abs(r), 2) +
np.power(np.abs(s), 2)))
elif self.dim == 3:
dx, dy, dz = op.grad(sol, dim=self.dim, **self.kwargs)
r -= 1. / (12. * T * mt**2) * dx
s -= 1. / (12. * T * mt**2) * dy
k -= 1. / (12. * T * mt**2) * dz
weights = np.maximum(1, np.sqrt(np.power(np.abs(r), 2) +
np.power(np.abs(s), 2) +
np.power(np.abs(k), 2)))
elif self.dim == 4:
dx, dy, dz, dt = op.grad(sol, dim=self.dim, **self.kwargs)
r -= 1. / (16 * T * mt**2) * dx
s -= 1. / (16 * T * mt**2) * dy
k -= 1. / (16 * T * mt**2) * dz
u -= 1. / (16 * T * mt**2) * dt
weights = np.maximum(1, np.sqrt(np.power(np.abs(r), 2) +
np.power(np.abs(s), 2) +
np.power(np.abs(k), 2) +
np.power(np.abs(u), 2)))
# FISTA update
t = (1 + np.sqrt(4 * told**2)) / 2.
if self.dim >= 1:
p = r / weights
r = p + (told - 1) / t * (p - pold)
pold = p
rr = deepcopy(r)
if self.dim >= 2:
q = s / weights
s = q + (told - 1) / t * (q - qold)
ss = deepcopy(s)
qold = q
if self.dim >= 3:
o = k / weights
k = o + (told - 1) / t * (o - kold)
kk = deepcopy(k)
kold = o
if self.dim >= 4:
m = u / weights
u = m + (told - 1) / t * (m - uold)
uu = deepcopy(u)
uold = m
told = t
iter += 1
try:
type(crit) == str
except NameError:
crit = "MAX_IT"
t_end = time()
exec_time = t_end - t_init
if self.verbosity in ['HIGH', 'ALL']:
print("Prox_TV: obj = {0}, rel_obj = {1}, {2}, iter = {3}".format(
obj, rel_obj, crit, iter))
print("exec_time = ", exec_time)
return sol
class proj(func):
r"""
Base class which defines the attributes of the `proj` objects.
See generic attributes descriptions of the
:class:`pyunlocbox.functions.func` base class.
Parameters
----------
epsilon : float, 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'.
Notes
-----
* All indicator functions (projections) evaluate to zero by definition.
"""
def __init__(self, epsilon=1, method='FISTA', **kwargs):
super(proj, self).__init__(**kwargs)
self.epsilon = epsilon
self.method = method
def _eval(self, x):
# Matlab version returns a small delta to avoid division by 0 when
# evaluating relative tolerance. Here the delta is added in the solve
# function if the sum of the objective functions is zero.
# np.spacing(1.0) is equivalent to matlab eps = eps(1.0)
# return np.spacing(1.0)
return 0
class proj_b2(proj):
r"""
L2-ball function object.
This function is the indicator function :math:`i_S(z)` of the set S which
is zero if `z` is in the set and infinite otherwise. The set S is defined
by :math:`\left\{z \in \mathbb{R}^N \mid \|A(z)-y\|_2 \leq \epsilon
\right\}`.
See generic attributes descriptions of the
:class:`pyunlocbox.functions.proj` base class. Note that the constructor
takes keyword-only parameters.
Notes
-----
* The `tol` parameter is defined as the tolerance for the projection on the
L2-ball. The algorithm stops if :math:`\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
:math:`\operatorname{arg\,min}\limits_z \frac{1}{2} \|x-z\|_2^2 + i_S(z)`
which has an identical solution as
:math:`\operatorname{arg\,min}\limits_z \|x-z\|_2^2` such that
:math:`\|A(z)-y\|_2 \leq \epsilon`. It is thus a projection of the vector
`x` onto an L2-ball of diameter `epsilon`.
Examples
--------
>>> import pyunlocbox
>>> f = pyunlocbox.functions.proj_b2(y=[1, 1])
>>> x = [3, 3]
>>> f.eval(x)
0
>>> f.prox(x, 0)
array([ 1.70710678, 1.70710678])
"""
def __init__(self, **kwargs):
# Constructor takes keyword-only parameters to prevent user errors.
super(proj_b2, self).__init__(**kwargs)
def _prox(self, x, T):
crit = None # Stopping criterion.
niter = 0 # Number of iterations.
# Tight frame.
if self.tight:
tmp1 = self.A(x) - self.y()
scale = self.epsilon / np.sqrt(np.sum(tmp1 * tmp1, axis=0))
tmp2 = tmp1 * np.minimum(1, scale) # Scaling.
sol = x + self.At(tmp2 - tmp1) / self.nu
crit = 'TOL'
u = np.nan
# Non tight frame.
else:
# Initialization.
sol = x
u = np.zeros(np.shape(self.y()))
if self.method is 'FISTA':
v_last = u
t_last = 1.
elif self.method is not 'ISTA':
raise ValueError('The method should be either FISTA or ISTA.')
# Tolerance around the L2-ball.
epsilon_low = self.epsilon / (1. + self.tol)
epsilon_up = self.epsilon / (1. - self.tol)
# Check if we are already in the L2-ball.
norm_res = np.linalg.norm(self.y() - self.A(sol), 2)
if norm_res <= epsilon_up:
crit = 'INBALL'
# Projection onto the L2-ball
while not crit:
niter += 1
# Residual.
res = self.A(sol) - self.y()
norm_res = np.linalg.norm(res, 2)
if self.verbosity is 'HIGH':
print(' proj_b2 iteration {:3d}: epsilon = {:.2e}, '
'||y-A(z)||_2 = {:.2e}'.format(niter, self.epsilon,
norm_res))
# Scaling for projection.
res += u * self.nu
norm_proj = np.linalg.norm(res, 2)
ratio = min(1, self.epsilon / norm_proj)
v = 1. / self.nu * (res - res * ratio)
if self.method is 'FISTA':
t = (1. + np.sqrt(1. + 4. * t_last**2.)) / 2. # Time step.
u = v + (t_last - 1.) / t * (v - v_last)
v_last = v
t_last = t
else:
u = v
# Current estimation.
sol = x - self.At(u)
# Stopping criterion.
if norm_res >= epsilon_low and norm_res <= epsilon_up:
crit = 'TOL'
elif niter >= self.maxit:
crit = 'MAXIT'
if self.verbosity in ['LOW', 'HIGH']:
norm_res = np.linalg.norm(self.y() - self.A(sol), 2)
print(' proj_b2: epsilon = {:.2e}, ||y-A(z)||_2 = {:.2e}, '
'{}, niter = {}'.format(self.epsilon, norm_res, crit,
niter))
return sol