/
KSD.py
699 lines (543 loc) · 23.4 KB
/
KSD.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
#!/usr/bin/env python3
# -*- coding: utf-8 -*-
"""
Created on Wed Aug 8 10:41:21 2018
@author: Anna SONG
KSD.py implements the `2D Kendall Shape Dictionary' that we propose in
[Song, Uhlmann, Fageot & Unser, Dictionary Learning for 2D Kendall Shapes (2019)].
Given a dataset of preshapes z_1,...,z_K in (C^N,Phi) endowed with the Hermitian
product Phi, the original problem to be minimized by
the dictionary D = [d_1,...,d_J] \in \C^{N \times J}
the weights A = [a_1,...,a_K] \in \C^{J \times K}
is written as
min_{D,A} sum_k d_F([z_k],[D a_k])^2 (*)
- d_F is the Full Procrustes distance
- where the atoms d_j in C^N are constrained to be preshapes
- the weights a_k in C^J are subject to |a_k|_0 <= N_0
- and such that D a_k are preshapes.
As shown in our article, this can be reformulated as as simple L2 problem in the complex setting
min_{D,A} sum_k || z_k - D a_k ||_Phi^2 (**)
- where the atoms d_j in C^N are constrained to be preshapes
- and the weights a_k in C^J are subject to |a_k|_0 <= N_0.
The solutions to the two minimization problems are equivalent up to a rescaling of the weights a_k.
KSD.py implements the minimization in the form of (**).
The 2D Kendall Shape Dictionary classically alternates between:
- a sparse coding step : the weights A are updated using a Cholesky-based
Order Recursive Matching Pursuit (ORMP), as a direct adaptation to the
complex setting of Mairal's implementation for the real setting in the SPAMS toolbox.
- a dictionary update : following the Method of Optimal Directions (MOD),
we update D as
D <- [z_1,...,z_K] @ A^H @ (A @ A^H)^{-1}
D <- Pi_S(D) (center and normalize all the atoms d_j)
Rmks:
- there is no need for a gradient descent with respect to D, although the expression
of the gradient can be explicited:
Grad_D [ sum_k || z_k - D a_k ||_Phi^2 ]
= 2 Phi (D @ sum_k a_k @ a_k^H - sum_k z_k @ a_k^H )
- after the sparse coding step in (**), D a_k are not preshapes, but just centered.
One can decide to rescale them so as to obtain preshapes.
- the sparse coding step ensures that the reconstructions D a_k in (**) are optimally
aligned to the data z_k (resp., optimally rotated in (*) if considering the corresponding
preshapes).
- in the dictionary update step:
- if A @ A^H is not invertible, we use the
SVD decomposition of A (see the article).
- In practice, the preshaping step D <- Pi_S(D) is just the normalization of
all the columns d_j (because the previous operation already centers them,
since z_1,...,z_K are already centered).
"""
import numpy as np
from numpy import pi
import matplotlib.pyplot as plt
import os
import scipy
import random
import time
from settings import M, N, K, data_type, test_k_set, database, choice, shapes # dataset contained in shapes
from settings import draw, drawMany
from settings import Psi, Phi, norm, her
from settings import multi_complex2real
from shapespace import theta, align_rot
import multiprocessing as mp
from itertools import product
from itertools import repeat
''' INITIALIZATION '''
def initializeD_c(J,dataset) :
'''Initialize the dictionary: randomly copy the existing data.'''
D_c = np.zeros((N,J),dtype = complex)
indices = np.arange(K) ; random.shuffle(indices)
for j in range(J) :
k = indices[j]
D_c[:,j] = dataset[k]
return D_c
''' SPARSE CODING '''
def ORMP_cholesky(D_c,dataset,N0) : # columns of D_c must be preshapes
'''Order Recursive Matching Pursuit with Cholesky-based optimisation,
as in the SPAMS toolbox of Mairal et al. (2009).
This is a direct adaptation of their code in C++ to the complex setting
with a Hermitian product Phi.
The columns of D_c must all be preshapes.
'''
# start = time.time()
J = D_c.shape[1]
K = len(dataset)
A_c = np.zeros((J,K),dtype = complex)
G = D_c.T.conj() @ Phi @ D_c
for k in range(K) :
z = dataset[k]
vM = np.zeros((N0),dtype = complex)
rM = np.zeros((N0),dtype=int)
Un = np.zeros((N0,N0), dtype = complex)
Undn = np.zeros((J,N0),dtype = complex)
Gs = np.zeros((N0,J),dtype = complex)
norm2 = np.ones(J)
rM[:] = -1
Rdn = D_c.T.conj() @ Phi @ z
scores = Rdn.copy()
for l in range(N0) :
currind = np.argmax(np.abs(scores))
if np.abs(scores[currind]) < 1e-8 :
#print('finished with a sparsity l =',l,'for data',k)
break
RU = scores[currind]
vM[l] = RU
rM[l] = currind
Un[l,l] = -1
Un[:l,l] = Un[:l,:l] @ Undn[currind,:l].T
Un[:,l] = - Un[:,l] / np.sqrt(norm2[currind])
if l == N0-1 :
break
Gs[l] = G[currind]
Undn[:,l] = Un[:l+1,l].T.conj() @ Gs[:l+1,:]
Rdn = Rdn - vM[l] * Undn[:,l].conj()
norm2 = norm2 - np.abs(Undn[:,l])**2
norm2[norm2<0] = 0 # sometimes happens to have small negative numbers
non_null = (norm2 > 1e-8)
scores = np.zeros(J,dtype=complex)
scores[non_null] = Rdn[non_null] / np.sqrt(norm2)[non_null]
vM = Un @ vM
indices = [j for j in rM if j >= 0]
A_c[indices,k] = vM[:len(indices)]
# elapsed = (time.time() - start)
# print('time : ',elapsed)
return A_c
'''ORMP_cholesky_real() and OMP() are not used in this script, but it is interesting to
write them as a comparison to ORMP_cholesky(). Empirically, OMP() behaves badly compared
ORMP_cholesky_real() for small values of N0, and becomes better after increasing N0.
ORMP_cholesky_real() is less performing than ORMP_cholesky() because the former
can only rely on real weights in sparse coding, contrarily to the latter. '''
# NOT USED HERE
def ORMP_cholesky_real(DM,N0,J,dataset) :
''' Same as ORMP_cholesky but in the real setting.
It copies the function of Mairal et al. (2009) for the SPAMS toolbox, as a comparison. '''
sqrtPsi = scipy.linalg.sqrtm(Psi).real
import spams
K = len(dataset)
dataset_r = multi_complex2real(dataset)
DR = np.asfortranarray(sqrtPsi @ multi_complex2real(DM.T).T)
X = np.asfortranarray(sqrtPsi @ dataset_r.T)
AM = spams.omp(X,D=DR,L=N0)
AM = np.array(AM.todense()).reshape(J,K)
return AM
# NOT USED HERE
def matPlus(B) :
'''Auxiliary function of OMP().
Gives the pseudo-inverse (B^H @ Phi @ B)^{-1} @ B^H @ Phi of a matrix B
with respect to the Hermitian product Phi.'''
try :
Bp = np.linalg.inv(B.T.conj() @ Phi @ B) @ B.T.conj() @ Phi
except np.linalg.LinAlgError :
global prob
prob = B.T.conj() @ Phi @ B
print('not invertible !!!\n',prob)
else :
return Bp
# NOT USED HERE
def OMP(D_c,dataset,N0) : # columns of D_c must be preshapes
''' Complex Orthogonal Matching Pursuit (complex OMP version).
Requires that the atoms contained in D_c are all preshapes.
'''
print('Caution! this is a naive implementation that does not handle singular cases.')
J = D_c.shape[1]
K = len(dataset)
A_c = np.zeros((J,K),dtype = complex)
for k in range(K) :
z = dataset[k]
rem = z
indices = np.array([],dtype = int)
for i in range(N0) :
corr = D_c.T.conj() @ Phi @ rem
corr[indices] = 0 # mathematically should be zero
j0 = int(np.argmax(np.abs(corr)))
indices = np.append(indices,j0)
DD_I = D_c[:,indices]
DD_I_plus = matPlus(DD_I)
P = DD_I @ DD_I_plus
rem = z - P @ z
a_I = DD_I_plus @ z # in C^N0
A_c[indices,k] = a_I
return A_c
''' PRESHAPING THE ATOMS OF D_c '''
# NOT USED HERE
def proj_S(D_c) : # columns of D_c are general configurations in C^N
'''converts non-zero columns of D_c to preshapes, by projecting them on the subspace
rthogonal to u {d | u.T.conj() Phi d = 0)}
and then normalizing them. Those which are zero are left as they are.
'''
from settings import uu
J = D_c.shape[1]
for j in range(J) :
d = D_c[:,j]
no = norm(d)
if no > 1e-8 :
orth_d = her(uu,d)*uu
d = d - orth_d
d = d / norm(d)
D_c[:,j] = d
return D_c
def normalize(D_c) : # columns of D_c are centered configurations
''' Same as proj_S(), but we know that all the columns are already centered. '''
J = D_c.shape[1]
for j in range(J) :
d = D_c[:,j]
no = norm(d)
if no > 1e-8 :
D_c[:,j] = d/no
return D_c
''' DISPLAY FUNCTIONS '''
from settings import z0
def display_res(dataset,DA,k,save = False,directory = None) :
'''Auxiliary function of display. display_res() show the reconstruction for data z_k'''
print('k = ',k)
z = dataset[k]
w = DA[:,k] # not necessarily normalized
ta = theta(z,z0) ; z = np.exp(1j*ta)*z ; w = np.exp(1j*ta)*w
norm_error = np.round(norm(z - w).real,4)
draw(z,i = 7,comparing = True) # i = 7 for gray
draw(w,i = 0)
plt.axis('equal')
plt.axis('off')
title_text = '$ |z_k - D \\alpha_k|_\\Phi : {} \% $'.format(np.round(100*norm_error,4))
plt.title(title_text, fontsize = 20, y = 1)
if save :
plt.savefig(directory + '/rec_'+str(k)+'.png',dpi = 200)
plt.show()
def display(D_c,A_c,dataset,test_k_set=test_k_set,save = False,directory = None) :
'''Displays some reconstructions, and optionally the sparse
coefficients A.'''
DA = D_c @ A_c
for k in test_k_set :
display_res(dataset,DA,k,save = save,directory = directory)
''' MAIN ALGORITHM : 2D Kendall Shape Dictionary with the Method of Optimal Directions and
Cholesky-optimized Order Recursive Matching Pursuit '''
def reciprocal(sigmas):
'''Auxiliary function of KSD_optimal_directions().
Given a 1D array of non-negative elements called sigmas, with possibly zero elements,
returns the array of multiplicative inverses whenever possible, and leaves the zeroes.'''
sigmas_rec = np.zeros_like(sigmas)
for i,x in enumerate(sigmas) :
if x != 0 :
sigmas_rec[i] = 1/x
return sigmas_rec
def fill_diagonal(sigmas_rec,J,K):
'''Auxiliary function of KSD_optimal_directions().
Fills in the diagonal of a matrix of shape (K,J).'''
Sigma_rec = np.zeros((K,J))
np.fill_diagonal(Sigma_rec,sigmas_rec)
return Sigma_rec
def KSD_optimal_directions(dataset,N0,J,init=None,Ntimes = 100,verbose=False,save=False,directory = None) :
''' The 2D Kendall Shape Dictionary classically alternates between:
- a sparse coding step : the weights A are updated using a Cholesky-based
Order Recursive Matching Pursuit (ORMP), as a direct adaptation to the
complex setting of Mairal's implementation for the real setting in the SPAMS toolbox.
- a dictionary update : following the Method of Optimal Directions (MOD),
we update D as
D <- [z_1,...,z_K] @ A^H @ (A @ A^H)^{-1}
D <- Pi_S(D) (center and normalize all the non-null atoms d_j)
and then replace under-utilized or null atoms by randomly picked data.
An atom d_j is arbitrarily said to be under-utilized if
(nb of data using d_j) / (K*N0) < 1 / (50*J)
Parameters:
- dataset in C^{(K,n)} is a complex array containing the horizontally stacked dataset [z_1,...,z_K]^T
- N0 determines the L0 sparsity of the weights a_k
- J fixes the number of atoms that we want to learn
- init = None initializes the dictionary with randomly picked data shapes.
if init is a given (n,J) complex array, then the initialization starts with init.
- Ntimes is the number of iterations
- if verbose == True, the algorithm keeps track of the loss function E to be minimized at each iteration.
It saves time to set verbose = False.
'''
K = len(dataset)
if type(init) == np.ndarray :
D_c = init
else :
D_c = initializeD_c(J,dataset)
if verbose :
print('Initializing the dictionary.')
drawMany(D_c.T,force = True)
lossCurve = np.array([])
print("Let's wait for {} iterations...".format(Ntimes))
start = time.time()
for t in range(Ntimes) :
if t % 5 == 0 :
print('t =',t)
A_c = ORMP_cholesky(D_c,dataset,N0)
if verbose :
diffs = dataset.T - D_c @ A_c
if K < 10000 :
E = np.diag(diffs.T.conj() @ Phi @ diffs).sum().real
else :
E = 0
for k in range(K) :
E += (diffs[:,k].conj() @ Phi @ diffs[:,k]).real
lossCurve = np.append(lossCurve,E)
try :
Mat = np.linalg.inv(A_c @ A_c.T.conj())
except np.linalg.LinAlgError :
global A_error
A_error = A_c
print('A @ A^H not invertible, using SVD')
U,sigmas,VH = np.linalg.svd(A_c)
sigmas_rec = reciprocal(sigmas)
Sigma_rec = fill_diagonal(sigmas_rec,J,K)
D_c = dataset.T @ VH.T.conj() @ Sigma_rec @ U.T.conj()
else :
D_c = dataset.T @ A_c.T.conj() @ Mat
D_c = normalize(D_c) # the new atoms are preshaped
purge_j = np.where((np.abs(A_c)>1e-3).sum(axis=1)/K < N0/(5*J))[0]
for j in range(J) :
if norm(D_c[:,j]) < 1e-8 or j in purge_j :
print('purged ',j,'at iteration',t)
D_c[:,j] = shapes[np.random.randint(K)]
print('computing the final weights...')
A_c = ORMP_cholesky(D_c,dataset,N0) ;
elapsed = (time.time() - start)
print('duration of the algorithm: ', np.round(elapsed,2), 'seconds')
diffs = dataset.T - D_c @ A_c
if K < 10000 :
E = np.diag(diffs.T.conj() @ Phi @ diffs).sum().real
else :
E = 0
for k in range(K) :
E += (diffs[:,k].conj() @ Phi @ diffs[:,k]).real
print('FINAL RESULTS')
if verbose :
display(D_c,A_c,dataset, save = save,directory=directory)
lossCurve = np.append(lossCurve,E)
plt.figure() ;
plt.plot(np.arange(len(lossCurve)),lossCurve)
plt.title('Loss curve for the KSD algorithm')
if save : plt.savefig(directory + '/losscurve.png',dpi = 100)
plt.show()
drawMany(D_c.T,force = True,show = False)
plt.title('KSD dictionary N0 = {} J = {}'.format(N0,J))
if save :
plt.savefig(directory + '/dico_KSD.png',dpi = 200)
plt.show()
D_al = align_rot(D_c.T).T
drawMany(D_al.T,force = True,show = False)
plt.title('KSD rotated dictionary N0 = {} J = {}'.format(N0,J))
if save :
plt.savefig(directory + '/dico_KSD_rotated.png',dpi = 200)
plt.show()
print('Final loss : ',E)
print('RMSE :',np.sqrt(E/K))
if save :
text_file = open(directory + '/readme.txt','a')
text_file.write('\nduration of the algorithm: {} s \n \n'.format(np.round(elapsed,2)))
text_file.write('Final loss: {}\n'.format(E))
text_file.close()
return D_c,A_c,E
def OMRP_multiproc_helper(k,D_c,dataset,G,N0) :
J = D_c.shape[1]
z = dataset[k]
vM = np.zeros((N0),dtype = complex)
rM = np.zeros((N0),dtype=int)
Un = np.zeros((N0,N0), dtype = complex)
Undn = np.zeros((J,N0),dtype = complex)
Gs = np.zeros((N0,J),dtype = complex)
norm2 = np.ones(J)
rM[:] = -1
Rdn = D_c.T.conj() @ Phi @ z
scores = Rdn.copy()
for l in range(N0) :
currind = np.argmax(np.abs(scores))
if np.abs(scores[currind]) < 1e-8 :
#print('finished with a sparsity l =',l,'for data',k)
break
RU = scores[currind]
vM[l] = RU
rM[l] = currind
Un[l,l] = -1
Un[:l,l] = Un[:l,:l] @ Undn[currind,:l].T
Un[:,l] = - Un[:,l] / np.sqrt(norm2[currind])
if l == N0-1 :
break
Gs[l] = G[currind]
Undn[:,l] = Un[:l+1,l].T.conj() @ Gs[:l+1,:]
Rdn = Rdn - vM[l] * Undn[:,l].conj()
norm2 = norm2 - np.abs(Undn[:,l])**2
norm2[norm2<0] = 0 # sometimes happens to have small negative numbers
non_null = (norm2 > 1e-8)
scores = np.zeros(J,dtype=complex)
scores[non_null] = Rdn[non_null] / np.sqrt(norm2)[non_null]
vM = Un @ vM
indices = [j for j in rM if j >= 0]
a = np.zeros(J,dtype = complex)
a[indices] = vM[:len(indices)]
return a
def KSD_optimal_directions_multiproc_ORMP(dataset,N0,J,init=None,Ntimes = 100,batch_size=1024,verbose=False,save=False,directory = None) :
'''See KSD_optimal_directions().
THIS FUNCTION IS INTERESTING ONLY FOR LARGE DATASETS (K > 4000 for instance).
In this function, the ORMP sparse coding step in computed in parallel and independently
on the different z_k, by randomly chosen batches of size given in batch_size,
thanks to the multiprocessing library run with the function OMRP_multiproc_helper().
'''
K = len(dataset)
if type(init) == np.ndarray :
D_c = init
else :
D_c = initializeD_c(J,dataset)
if verbose :
print('Initializing the dictionary.')
drawMany(D_c.T,force = True)
lossCurve = np.array([])
A_c = np.zeros((J,K),dtype = complex)
print("Let's wait for {} iterations...".format(Ntimes))
start = time.time()
for t in range(Ntimes) :
if t % 5 == 0 :
print('t =',t)
indices = np.arange(K)
random.shuffle(indices)
indices = indices[:batch_size]
'parallel ORMP used, avoids distorted atoms but increases run-time'
G = D_c.T.conj() @ Phi @ D_c
pool = mp.Pool(mp.cpu_count())
results = pool.starmap(OMRP_multiproc_helper,zip(indices,repeat(D_c),repeat(dataset),repeat(G),repeat(N0)))
pool.close()
A_c[:,indices] = np.array(results).T
if verbose :
diffs = dataset.T - D_c @ A_c
if K < 10000 :
E = np.diag(diffs.T.conj() @ Phi @ diffs).sum().real
else :
E = 0
for k in range(K) :
E += (diffs[:,k].conj() @ Phi @ diffs[:,k]).real
lossCurve = np.append(lossCurve,E)
try :
Mat = np.linalg.inv(A_c @ A_c.T.conj())
except np.linalg.LinAlgError :
global A_error
A_error = A_c
print('A @ A^H not invertible, using SVD')
U,sigmas,VH = np.linalg.svd(A_c)
sigmas_rec = reciprocal(sigmas)
Sigma_rec = fill_diagonal(sigmas_rec,J,K)
D_c = dataset.T @ VH.T.conj() @ Sigma_rec @ U.T.conj()
else :
D_c = dataset.T @ A_c.T.conj() @ Mat
D_c = normalize(D_c) # the new atoms are preshaped
purge_j = np.where((np.abs(A_c)>1e-3).sum(axis=1)/K < N0/(5*J))[0]
purged_list = []
for j in range(J) :
if norm(D_c[:,j]) < 1e-8 or j in purge_j :
purged_list += [j]
#print('purged ',j,'at iteration',t)
D_c[:,j] = shapes[np.random.randint(K)]
if len(purged_list) > 0 :
print('purged atoms ', purged_list, 'at iteration',t)
print('using parallel ORMP to compute the final weights...')
'parallel ORMP used for the final computation'
G = D_c.T.conj() @ Phi @ D_c
pool = mp.Pool(mp.cpu_count())
results = pool.starmap(OMRP_multiproc_helper,zip(range(K),repeat(D_c),repeat(dataset),repeat(G),repeat(N0)))
pool.close()
A_c = np.array(results).T
elapsed = (time.time() - start)
print('duration of the algorithm: ', np.round(elapsed,2), 'seconds')
diffs = dataset.T - D_c @ A_c
if K < 10000 :
E = np.diag(diffs.T.conj() @ Phi @ diffs).sum().real
else :
E = 0
for k in range(K) :
E += (diffs[:,k].conj() @ Phi @ diffs[:,k]).real
print('FINAL RESULTS')
display(D_c,A_c,dataset, save = save,directory=directory)
if verbose :
lossCurve = np.append(lossCurve,E)
plt.figure() ;
plt.plot(np.arange(len(lossCurve)),lossCurve)
plt.title('Loss curve for the KSD algorithm')
if save : plt.savefig(directory + '/losscurve.png',dpi = 100)
plt.show()
drawMany(D_c.T,force = False,show = False)
plt.title('KSD dictionary N0 = {} J = {}'.format(N0,J))
if save :
plt.savefig(directory + '/dico_KSD.png',dpi = 200)
plt.show()
D_al = align_rot(D_c.T).T
drawMany(D_al.T,force = False,show = False)
plt.title('KSD rotated dictionary N0 = {} J = {}'.format(N0,J))
if save :
plt.savefig(directory + '/dico_KSD_rotated.png',dpi = 200)
plt.show()
print('Final loss : ',E)
print('RMSE :',np.sqrt(E/K))
if save :
text_file = open(directory + '/readme.txt','a')
text_file.write('\nduration of the algorithm: {} s \n \n'.format(np.round(elapsed,2)))
text_file.write('Final loss: {}\n'.format(E))
text_file.write('Final RMSE: {}\n'.format(np.sqrt(E/K)))
text_file.close()
return D_c,A_c,E
def KSD_RMSE_curve(dataset,J) :
'''
For a fixed J (nb of atoms imposed), for any N0 <= J,
run 2DKSD and record the RMSE(N0,J). Return the RMSE curve.
Do not run on too large datasets!! (K > 10000)
'''
K = len(dataset)
curve = np.zeros(J)
for N0 in range(1,J+1) :
D_c,A_c,E_KSD = KSD_optimal_directions(shapes,N0,J,
Ntimes = 30,verbose=False)
curve[N0-1] = np.sqrt(E_KSD/K)
plt.plot(curve)
plt.show()
return curve
if __name__ == '__main__' :
learn_dico = True
N0,J = 5,10
SAVE = True
MULTIPROC_ORMP = False
# True if you want to multiprocess ORMP
# interesting only for large datasets K > 4000
if learn_dico :
directory = 'RESULTS/' + database[choice] + '/N0_'+ str(N0) + '_J_' + str(J)
if SAVE :
if not os.path.exists(directory):
print('CREATING THE FOLDER')
os.makedirs(directory)
else :
if os.path.exists(directory + '/dico_KSD.png') :
SAVE = False
print('WILL NOT OVERWRITE PREVIOUS SAVED RESULTS')
if MULTIPROC_ORMP :
from functools import partial
KSD = partial(KSD_optimal_directions_multiproc_ORMP,
batch_size = K) # 1024 ~ 2048 ~ 4096
# you can set the batch size to the size of the dataset to
#multiprocess ORMP on the whole dataset
print('USING PARALLELIZED AND/OR BATCHED VERSION OF ORMP')
print('change the batch size if needed')
else :
KSD = KSD_optimal_directions
D_c,A_c,E_KSD = KSD(shapes,N0,J,init = None,
Ntimes = 30,verbose=False,
save=SAVE,directory=directory)
# change Ntimes as needed if you see that
# the loss curve has not finished converging (in verbose=True mode)
# and put verbose = False to speed up the process and have actual run-time