-
Notifications
You must be signed in to change notification settings - Fork 21
/
cp_apr.m
1691 lines (1440 loc) · 60.2 KB
/
cp_apr.m
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
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
function [M, Minit, output] = cp_apr(X, R, varargin)
%CP_APR Compute nonnegative CP with alternating Poisson regression.
%
% M = CP_APR(X, R) computes an estimate of the best rank-R CP model of a
% nonnegative tensor X using an alternating Poisson regression. This is
% most appropriate for sparse count data (i.e., nonnegative integer
% values) because it uses Kullback-Liebler divergence. The input X can
% be a tensor or sptensor. The result M is a ktensor. Input data must be
% nonnegative, and the computed ktensor factors are all nonnegative.
%
% Different algorithm variants are available (selected by the 'alg'
% parameter):
% 'pqnr' - row subproblems by projected quasi-Newton (default)
% 'pdnr' - row subproblems by projected damped Hessian
% 'mu' - multiplicative update (default in version 2.5)
%
% M = CP_APR(X, R, 'param', value, ...) specifies optional parameters and
% values. Some parameters work in all situations, others apply only for
% a particular choice of algorithm.
%
% Valid parameters and their default values are:
% 'alg' - Algorithm ['mu'|'pdnr'|'pqnr'] {'pqnr'}
% 'stoptol' - Tolerance on the overall KKT violation {1.0e-4}
% 'stoptime' - Maximum number of seconds to run {1e6}
% 'maxiters' - Maximum number of iterations {1000}
% 'init' - Initial guess [{'random'}|ktensor]
% 'maxinneriters' - Maximum inner iterations per outer iteration {10}
% 'epsDivZero' - Safeguard against divide by zero {1.0e-10}
% 'printitn' - Print every n outer iterations; 0 for none {1}
% 'printinneritn' - Print every n inner iterations {0}
%
% Additional input parameters for algorithm 'mu':
% 'kappa' - Offset to fix complementary slackness {100}
% 'kappatol' - Tolerance on complementary slackness {1.0e-10}
%
% Additional input parameters for algorithm 'pdnr':
% 'epsActive' - Bertsekas tolerance for active set {1.0e-8}
% 'mu0' - Initial damping parameter {1.0e-5}
% 'precompinds' - Precompute sparse tensor indices {true}
% 'inexact' - Compute inexact Newton steps {true}
%
% Additional input parameters for algorithm 'pqnr':
% 'epsActive' - Bertsekas tolerance for active set {1.0e-8}
% 'lbfgsMem' - Number vector pairs to store for L-BFGS {3}
% 'precompinds' - Precompute sparse tensor indices {true}
%
% [M,M0] = CP_APR(...) also returns the initial guess.
%
% [M,M0,out] = CP_APR(...) also returns additional output:
% out.kktViolations - maximum KKT violation per iteration
% out.nInnerIters - number of inner iterations per outer iteration
% out.obj - final negative log-likelihood objective
% out.ttlTime - time algorithm took to converge or reach max
% out.times - cumulative time through each outer iteration
% If algorithm is 'mu':
% out.nViolations - number of factor matrices needing complementary
% slackness adjustment per iteration
% If algorithm is 'pdnr' or 'pqnr':
% out.nZeros - number of zero factor entries per iteration
%
% REFERENCES:
% * E. C. Chi and T. G. Kolda. On Tensors, Sparsity, and Nonnegative
% Factorizations, SIAM J. Matrix Analysis, 33(4):1272-1299, Dec. 2012,
% http://dx.doi.org/10.1137/110859063
% * S. Hansen, T. Plantenga and T. G. Kolda, Newton-Based Optimization
% for Kullback-Leibler Nonnegative Tensor Factorizations,
% Optimization Methods and Software, 2015,
% http://dx.doi.org/10.1080/10556788.2015.1009977
%
% See also CP_ALS, KTENSOR, TENSOR, SPTENSOR.
%
%MATLAB Tensor Toolbox.
%Copyright 2015, Sandia Corporation.
% This is the MATLAB Tensor Toolbox by T. Kolda, B. Bader, and others.
% http://www.sandia.gov/~tgkolda/TensorToolbox.
% Copyright (2015) Sandia Corporation. Under the terms of Contract
% DE-AC04-94AL85000, there is a non-exclusive license for use of this
% work by or on behalf of the U.S. Government. Export of this data may
% require a license from the United States Government.
% The full license terms can be found in the file LICENSE.txt
%% Set the algorithm choice and initial guess from input or defaults.
params = inputParser;
params.addParameter('alg', 'pqnr', @(x) (ismember(x,{'mu','pdnr','pqnr'})) );
params.addParameter('init','random', @(x) (isa(x,'ktensor') || ismember(x,{'random'})) );
params.KeepUnmatched = true;
params.parse(varargin{:});
alg = params.Results.alg;
Minit = params.Results.init;
% Extract the number of modes in tensor X.
N = ndims(X);
if (R <= 0)
error('Number of components requested must be positive');
end
%% Check that the data is nonnegative.
tmp = find(X < 0.0);
if (size(tmp,1) > 0)
error('Data tensor must be nonnegative for Poisson-based factorization');
end
%% Set up an initial guess for the factor matrices.
if isa(Minit,'ktensor')
% User provided an initial ktensor; validate it.
if (ndims(Minit) ~= N)
error('Initial guess does not have the right number of modes');
end
if (ncomponents(Minit) ~= R)
error('Initial guess does not have the right number of components');
end
for n = 1:N
if (size(Minit,n) ~= size(X,n))
error('Mode %d of the initial guess is the wrong size',n);
end
if (min(min(Minit.U{n})) < 0.0)
error('Initial guess has negative element in mode %d',n);
end
end
if (min(Minit.lambda) < 0.0)
error('Initial guess has a negative ktensor weight');
end
elseif strcmp(Minit,'random')
% Choose random values for each element in the range (0,1).
F = cell(N,1);
for n = 1:N
F{n} = rand(size(X,n),R);
end
Minit = ktensor(F);
end
%% Call a solver based on the choice of algorithm parameter, passing
% all the other input parameters.
if strcmp(alg,'mu')
[M, output] = tt_cp_apr_mu (X, R, Minit, params.Unmatched);
output.params.alg = 'mu';
elseif strcmp(alg,'pdnr')
[M, output] = tt_cp_apr_pdnr (X, R, Minit, params.Unmatched);
output.params.alg = 'pdnr';
elseif strcmp(alg,'pqnr')
[M, output] = tt_cp_apr_pqnr (X, R, Minit, params.Unmatched);
output.params.alg = 'pqnr';
end
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% Main algorithm PQNR
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function [M, out] = tt_cp_apr_pqnr(X, R, Minit, varargin)
%TT_CP_APR_PQNR Compute nonnegative CP with alternating Poisson regression.
%
% tt_cp_apr_pqnr(X, R, ...) computes an estimate of the best rank-R
% CP model of a tensor X using an alternating Poisson regression.
% The algorithm solves "row subproblems" in each alternating subproblem,
% using a quasi-Newton Hessian approximation.
% The function is typically called by cp_apr.
%
% The model is solved by nonlinear optimization, and the code literally
% minimizes the negative of log-likelihood. However, printouts to the
% console reverse the sign to show maximization of log-likelihood.
%
% The function call can specify optional parameters and values.
% Valid parameters and their default values are:
% 'stoptol' - Tolerance on the overall KKT violation {1.0e-4}
% 'stoptime' - Maximum number of seconds to run {1e6}
% 'maxiters' - Maximum number of iterations {1000}
% 'maxinneriters' - Maximum inner iterations per outer iteration {10}
% 'epsDivZero' - Safeguard against divide by zero {1.0e-10}
% 'printitn' - Print every n outer iterations; 0 for no printing {1}
% 'printinneritn' - Print every n inner iterations {0}
% 'epsActive' - Bertsekas tolerance for active set {1.0e-8}
% 'lbfgsMem' - Number vector pairs to store for L-BFGS {3}
% 'precompinds' - Precompute sparse tensor indices to run faster {true}
%
% Return values are:
% M - ktensor model with R components
% out.fnEvals - number of row obj fn evaluations per outer iteration
% out.kktViolations - maximum KKT violation per iteration
% out.nInnerIters - number of inner iterations per outer iteration
% out.nZeros - number of factor elements equal to zero per iteration
% out.obj - final log-likelihood objective
% (minimization objective is actually -1 times this)
% out.ttlTime - time algorithm took to converge or reach max
% out.times - cumulative time through each outer iteration
%
% REFERENCE: Samantha Hansen, Todd Plantenga, Tamara G. Kolda.
% Newton-Based Optimization for Nonnegative Tensor Factorizations,
% arXiv:1304.4964 [math.NA], April 2013,
% URL: http://arxiv.org/abs/1304.4964. Submitted for publication.
%
% See also CP_APR, KTENSOR, TENSOR, SPTENSOR.
%
%MATLAB Tensor Toolbox.
%Copyright 2015, Sandia Corporation.
% This is the MATLAB Tensor Toolbox by T. Kolda, B. Bader, and others.
% http://www.sandia.gov/~tgkolda/TensorToolbox.
% Copyright (2015) Sandia Corporation. Under the terms of Contract
% DE-AC04-94AL85000, there is a non-exclusive license for use of this
% work by or on behalf of the U.S. Government. Export of this data may
% require a license from the United States Government.
% The full license terms can be found in the file LICENSE.txt
%% Set algorithm parameters from input or by using defaults.
params = inputParser;
params.addParamValue('epsActive', 1e-8, @isscalar);
params.addParamValue('epsDivZero',1e-10,@isscalar);
params.addParamValue('lbfgsMem',3,@isscalar);
params.addParamValue('maxinneriters',10,@isscalar);
params.addParamValue('maxiters',1000,@(x) isscalar(x) & x > 0);
params.addParamValue('precompinds',true,@(x) isa(x,'logical'));
params.addParamValue('printinneritn',0,@isscalar);
params.addParamValue('printitn',1,@isscalar);
params.addParamValue('stoptime',1e6,@isscalar);
params.addParamValue('stoptol',1e-4,@isscalar);
params.parse(varargin{:});
%% Copy from params object.
epsActSet = params.Results.epsActive;
epsDivZero = params.Results.epsDivZero;
nSizeLBFGS = params.Results.lbfgsMem;
maxInnerIters = params.Results.maxinneriters;
maxOuterIters = params.Results.maxiters;
precomputeSparseIndices = params.Results.precompinds;
printInnerItn = params.Results.printinneritn;
printOuterItn = params.Results.printitn;
stoptime = params.Results.stoptime;
stoptol = params.Results.stoptol;
% Extract the number of modes in tensor X.
N = ndims(X);
% If the initial guess has any rows of all zero elements, then modify
% so the row subproblem is not taking log(0). Values will be restored to
% zero later if the unfolded X for the row has no nonzeros.
for n = 1:N
rowsum = sum(Minit{n},2);
tmpIx = find(rowsum == 0);
if (isempty(tmpIx) == false)
Minit{n}(tmpIx,1) = 1.0e-8;
end
end
% Start with the initial guess, normalized using the vector L1 norm.
M = normalize(Minit,[],1);
% Sparse tensor flag affects how Pi and Phi are computed.
if isa(X,'sptensor')
isSparse = true;
else
isSparse = false;
end
% Initialize output arrays.
fnEvals = zeros(maxOuterIters,1);
kktViolations = -ones(maxOuterIters,1);
nInnerIters = zeros(maxOuterIters,1);
nzeros = zeros(maxOuterIters,1);
times = zeros(maxOuterIters,1);
if (printOuterItn > 0)
fprintf('\nCP_PQNR (alternating Poisson regression using quasi-Newton)\n');
end
dispLineWarn = (printInnerItn > 0);
% Start the wall clock timer.
tic;
if (isSparse && precomputeSparseIndices)
% Precompute sparse index sets for all the row subproblems.
% Takes more memory but can cut execution time significantly in some cases.
if (printOuterItn > 0)
fprintf(' Precomputing sparse index sets...');
end
sparseIx = cell(N);
for n = 1:N
num_rows = size(M{n},1);
sparseIx{n} = cell(num_rows,1);
for jj = 1:num_rows
sparseIx{n}{jj} = find(X.subs(:,n) == jj);
end
end
if (printOuterItn > 0)
fprintf('done\n');
end
end
%% Main Loop: Iterate until convergence or a max threshold is reached.
for iter = 1:maxOuterIters
isConverged = true;
kktModeViolations = zeros(N,1);
countInnerIters = zeros(1,N);
% Alternate thru each factor matrix, A_1, A_2, ... , A_N.
for n = 1:N
% Shift the weight from lambda to mode n.
M = redistribute(M,n);
% Calculate Khatri-Rhao product of all matrices but the n-th.
if (isSparse == false)
% Data is not a sparse tensor.
Pi = tt_calcpi_prowsubprob(X, isSparse, M, R, n, N, []);
X_mat = double(tenmat(X,n));
end
num_rows = size(M{n},1);
isRowNOTconverged = zeros(1,num_rows);
% Loop over the row subproblems in mode n.
for jj = 1:num_rows
% Get data values for row jj of matricized mode n.
if (isSparse)
% Data is a sparse tensor.
if (precomputeSparseIndices == false)
sparse_indices = find(X.subs(:,n) == jj);
else
sparse_indices = sparseIx{n}{jj};
end
if (isempty(sparse_indices))
% The row jj of matricized tensor X in mode n is empty.
M{n}(jj,:) = 0;
continue
end
x_row = X.vals(sparse_indices);
% Calculate just the columns of Pi needed for this row.
Pi = tt_calcpi_prowsubprob(X, isSparse, M, ...
R, n, N, sparse_indices);
else
x_row = X_mat(jj,:);
end
% Get current values of the row subproblem variables.
m_row = M{n}(jj,:);
% Initialize L-BFGS storage for the row subproblem.
delm = zeros(R, nSizeLBFGS);
delg = zeros(R, nSizeLBFGS);
rho = zeros(nSizeLBFGS, 1);
lbfgsPos = 1;
m_rowOLD = [];
gradOLD = [];
% Iteratively solve the row subproblem with projected qNewton steps.
for i = 1:maxInnerIters
% Calculate the gradient.
[gradM, phi_row] = calc_grad(isSparse, Pi, epsDivZero, ...
x_row, m_row);
if (i == 1)
% Original cp_aprPQN_row code (and plb_row) does a gradient
% step to prime the L-BFGS approximation. However, it means
% a row subproblem that already converged wastes time
% doing a gradient step before checking KKT conditions.
% TODO: fix in a future release.
m_rowOLD = m_row;
gradOLD = gradM;
[m_row, f, f_unit, f_new, num_evals] ...
= tt_linesearch_prowsubprob(-gradM', gradM', ...
m_rowOLD, ...
1, 1/2, 10, 1.0e-4, ...
isSparse, x_row, Pi, ...
phi_row, dispLineWarn);
fnEvals(iter) = fnEvals(iter) + num_evals;
[gradM, phi_row] = calc_grad(isSparse, Pi, epsDivZero, ...
x_row, m_row);
end
% Compute the row subproblem kkt_violation.
% Experiments in the original paper used this:
%kkt_violation = norm(abs(min(m_row,gradM')),2);
% Now we use | KKT |_inf:
kkt_violation = max(abs(min(m_row,gradM')));
% Report largest row subproblem initial violation.
if ((i == 1) && (kkt_violation > kktModeViolations(n)))
kktModeViolations(n) = kkt_violation;
end
if (mod(i, printInnerItn) == 0)
fprintf(' Mode = %1d, Row = %d, InnerIt = %d', ...
n, jj, i);
if (i == 1)
fprintf(', RowKKT = %.2e\n', kkt_violation);
else
fprintf(', RowKKT = %.2e, RowObj = %.4e\n', ...
kkt_violation, -f_new);
end
end
% Check for row subproblem convergence.
if (kkt_violation < stoptol)
break;
else
% Not converged, so m_row will be modified.
isRowNOTconverged(jj) = 1;
end
% Update the L-BFGS approximation.
tmp_delm = m_row - m_rowOLD;
tmp_delg = gradM - gradOLD;
tmp_rho = 1 / (tmp_delm * tmp_delg);
if ((tmp_rho > 0.0) && (isinf(tmp_rho) == false))
delm(:,lbfgsPos) = tmp_delm;
delg(:,lbfgsPos) = tmp_delg;
rho(lbfgsPos) = tmp_rho;
else
% Rho is required to be positive; if not, then skip
% the L-BFGS update pair. The recommended safeguard for
% full BFGS is Powell damping, but not clear how to damp
% in 2-loop L-BFGS.
if (dispLineWarn)
fprintf('WARNING: skipping L-BFGS update, rho would be 1 / %.2e\n', ...
(tmp_delm * tmp_delg));
end
% Roll back lbfgsPos since it will increment later.
if (lbfgsPos == 1)
if (rho(nSizeLBFGS) > 0)
lbfgsPos = nSizeLBFGS;
else
% Fatal error, should not happen.
fprintf('ERROR: L-BFGS first iterate is bad\n');
return;
end
else
lbfgsPos = lbfgsPos - 1;
end
end
% Calculate the search direction.
search_dir = getSearchDirPqnr(m_row, gradM, epsActSet, ...
delm, delg, rho, lbfgsPos, ...
i, dispLineWarn);
lbfgsPos = mod(lbfgsPos, nSizeLBFGS) + 1;
m_rowOLD = m_row;
gradOLD = gradM;
% Perform a projected linesearch and update variables.
% Start from a unit step length, decrease by 1/2, stop with
% sufficient decrease of 1.0e-4 or at most 10 steps.
[m_row, f, f_unit, f_new, num_evals] ...
= tt_linesearch_prowsubprob(search_dir', gradOLD', m_rowOLD, ...
1, 1/2, 10, 1.0e-4, ...
isSparse, x_row, Pi, ...
phi_row, dispLineWarn);
fnEvals(iter) = fnEvals(iter) + num_evals;
end
M{n}(jj,:) = m_row;
countInnerIters(n) = countInnerIters(n) + i;
end
% Test if all row subproblems have converged, which means that
% no variables in this mode were changed.
if (sum(isRowNOTconverged) ~= 0)
isConverged = false;
end
% Shift weight from mode n back to lambda.
M = normalize(M,[],1,n);
% Total number of inner iterations for a given outer iteration,
% totalled across all modes and all row subproblems in each mode.
nInnerIters(iter) = nInnerIters(iter) + countInnerIters(n);
end
% Save output items for the outer iteration.
num_zero = 0;
for n = 1:N
num_zero = num_zero + nnz(find(M{n} == 0.0));
end
nzeros(iter) = num_zero;
kktViolations(iter) = max(kktModeViolations);
% Print outer iteration status.
if (mod(iter,printOuterItn) == 0)
fprintf('%4d. Ttl Inner Its: %d, KKT viol = %.2e, obj = %.8e, nz: %d\n', ...
iter, nInnerIters(iter), kktViolations(iter), tt_loglikelihood(X,M), ...
num_zero);
end
times(iter) = toc;
% Check for convergence
if (isConverged)
break;
end
if (times(iter) > stoptime)
fprintf('Exiting because time limit exceeded\n');
break;
end
end
t_stop = toc;
%% Clean up final result and set output items.
M = normalize(M,'sort',1);
loglike = tt_loglikelihood(X,M);
if (printOuterItn > 0)
% For legacy reasons, compute "fit", the fraction explained by the model.
% Fit is in the range [0,1], with 1 being the best fit.
normX = norm(X);
normresidual = sqrt( normX^2 + norm(M)^2 - 2 * innerprod(X,M) );
fit = 1 - (normresidual / normX);
fprintf('===========================================\n');
fprintf(' Final log-likelihood = %e \n', loglike);
fprintf(' Final least squares fit = %e \n', fit);
fprintf(' Final KKT violation = %7.7e\n', kktViolations(iter));
fprintf(' Total inner iterations = %d\n', sum(nInnerIters));
fprintf(' Total execution time = %.2f secs\n', t_stop);
end
out = struct;
out.params = params.Results;
out.obj = loglike;
out.kktViolations = kktViolations(1:iter);
out.fnEvals = fnEvals(1:iter);
out.nInnerIters = nInnerIters(1:iter);
out.nZeros = nzeros(1:iter);
out.times = times(1:iter);
out.ttlTime = t_stop;
end
%----------------------------------------------------------------------
function [grad_row, phi_row] = calc_grad(isSparse, Pi, eps_div_zero, x_row, m_row)
%function grad_row = calc_grad(isSparse, Pi, eps_div_zero, x_row, m_row)
% Compute the gradient for a PQNR row subproblem.
%
% isSparse - true if x_row is sparse, false if dense
% Pi - matrix
% eps_div_zero - safeguard value to prevent division by zero
% x_row - row vector of data values for the row subproblem
% m_row - vector of variables for the row subproblem
%
% Returns the gradient vector for a row subproblem.
if (isSparse)
v = m_row * Pi';
w = x_row' ./ max(v, eps_div_zero);
phi_row = w * Pi;
else
v = m_row * Pi';
w = x_row ./ max(v, eps_div_zero);
phi_row = w * Pi;
end
grad_row = (ones(size(phi_row)) - phi_row)';
end
%----------------------------------------------------------------------
function [d] = getSearchDirPqnr (m_row, grad, epsActSet, ...
delta_m, delta_g, rho, lbfgs_pos, ...
iters, disp_warn)
% Compute the search direction by projecting with L-BFGS.
%
% m_row - current variable values
% grad - gradient at m_row
% epsActSet - Bertsekas tolerance for active set determination
% delta_m - L-BFGS array of vector variable deltas
% delta_g - L-BFGS array of gradient deltas
% lbfgs_pos - pointer into L-BFGS arrays
%
% Returns
% d - search direction based on current L-BFGS and grad
%
% Adapted from MATLAB code of Dongmin Kim and Suvrit Sra written in 2008.
% Modified extensively to solve row subproblems and use a better linesearch;
% see the reference at the top of this file for details.
lbfgsSize = size(delta_m,2);
% Determine active and free variables.
% If epsActSet is zero, then the following works:
% fixedVars = find((m_row == 0) & (grad' > 0));
% For the general case this works but is less clear and assumes m_row > 0:
% fixedVars = find((grad' > 0) & (m_row <= min(epsActSet,grad')));
projGradStep = (m_row - grad') .* (m_row - grad' > 0);
wk = norm(m_row - projGradStep);
fixedVars = find((grad' > 0) & (m_row <= min(epsActSet,wk)));
d = -grad;
d(fixedVars) = 0;
if ((delta_m(:,lbfgs_pos)' * delta_g(:,lbfgs_pos)) == 0.0)
% Cannot proceed with this L-BFGS data; most likely the iteration
% has converged, so this is rarely seen.
if (disp_warn)
fprintf('WARNING: L-BFGS update is orthogonal, using gradient\n');
end
return;
end
alpha = ones(lbfgsSize,1);
k = lbfgs_pos;
% Perform an L-BFGS two-loop recursion to compute the search direction.
for i = 1 : min(iters, lbfgsSize)
alpha(k) = rho(k) * delta_m(:, k)' * d;
d = d - alpha(k) * delta_g(:, k);
k = lbfgsSize - mod(1 - k, lbfgsSize);
end
coef = 1 / rho(lbfgs_pos) / (delta_g(:, lbfgs_pos)' * delta_g(:, lbfgs_pos));
d = coef * d;
for i = 1 : min(iters, lbfgsSize)
k = mod(k, lbfgsSize) + 1;
b = rho(k) * delta_g(:, k)' * d;
d = d + (alpha(k) - b) * delta_m(:, k);
end
d(fixedVars) = 0;
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% Main algorithm PDNR
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function [M, out] = tt_cp_apr_pdnr(X, R, Minit, varargin)
%TT_CP_APR_PDNR Compute nonnegative CP with alternating Poisson regression.
%
% tt_cp_apr_pdnr(X, R, ...) computes an estimate of the best rank-R
% CP model of a tensor X using an alternating Poisson regression.
% The algorithm solves "row subproblems" in each alternating subproblem,
% using a Hessian of size R^2.
% The function is typically called by cp_apr.
%
% The model is solved by nonlinear optimization, and the code literally
% minimizes the negative of log-likelihood. However, printouts to the
% console reverse the sign to show maximization of log-likelihood.
%
% The function call can specify optional parameters and values.
% Valid parameters and their default values are:
% 'stoptol' - Tolerance on the overall KKT violation {1.0e-4}
% 'stoptime' - Maximum number of seconds to run {1e6}
% 'maxiters' - Maximum number of iterations {1000}
% 'maxinneriters' - Maximum inner iterations per outer iteration {10}
% 'epsDivZero' - Safeguard against divide by zero {1.0e-10}
% 'printitn' - Print every n outer iterations; 0 for no printing {1}
% 'printinneritn' - Print every n inner iterations {0}
% 'epsActive' - Bertsekas tolerance for active set {1.0e-8}
% 'mu0' - Initial damping parameter {1.0e-5}
% 'precompinds' - Precompute sparse tensor indices to run faster {true}
% 'inexact' - Compute inexact Newton steps {true}
%
% Return values are:
% M - ktensor model with R components
% out.fnEvals - number of row obj fn evaluations per outer iteration
% out.kktViolations - maximum KKT violation per iteration
% out.nInnerIters - number of inner iterations per outer iteration
% out.nZeros - number of factor elements equal to zero per iteration
% out.obj - final log-likelihood objective
% (minimization objective is actually -1 times this)
% out.ttlTime - time algorithm took to converge or reach max
% out.times - cumulative time through each outer iteration
%
% REFERENCE: Samantha Hansen, Todd Plantenga, Tamara G. Kolda.
% Newton-Based Optimization for Nonnegative Tensor Factorizations,
% arXiv:1304.4964 [math.NA], April 2013,
% URL: http://arxiv.org/abs/1304.4964. Submitted for publication.
%
% See also CP_APR, KTENSOR, TENSOR, SPTENSOR.
%
%MATLAB Tensor Toolbox.
%Copyright 2015, Sandia Corporation.
% This is the MATLAB Tensor Toolbox by T. Kolda, B. Bader, and others.
% http://www.sandia.gov/~tgkolda/TensorToolbox.
% Copyright (2015) Sandia Corporation. Under the terms of Contract
% DE-AC04-94AL85000, there is a non-exclusive license for use of this
% work by or on behalf of the U.S. Government. Export of this data may
% require a license from the United States Government.
% The full license terms can be found in the file LICENSE.txt
%% Set algorithm parameters from input or by using defaults.
params = inputParser;
params.addParamValue('epsActive', 1e-8, @isscalar);
params.addParamValue('epsDivZero',1e-10,@isscalar);
params.addParamValue('maxinneriters',10,@isscalar);
params.addParamValue('maxiters',1000,@(x) isscalar(x) & x > 0);
params.addParamValue('precompinds',true,@(x) isa(x,'logical'));
params.addParamValue('inexact',true,@(x) isa(x,'logical'));
params.addParamValue('mu0',1e-5,@isscalar);
params.addParamValue('printinneritn',0,@isscalar);
params.addParamValue('printitn',1,@isscalar);
params.addParamValue('stoptime',1e6,@isscalar);
params.addParamValue('stoptol',1e-4,@isscalar);
params.parse(varargin{:});
%% Copy from params object.
epsActSet = params.Results.epsActive;
epsDivZero = params.Results.epsDivZero;
maxInnerIters = params.Results.maxinneriters;
maxOuterIters = params.Results.maxiters;
mu0 = params.Results.mu0;
precomputeSparseIndices = params.Results.precompinds;
inexactNewton = params.Results.inexact;
printInnerItn = params.Results.printinneritn;
printOuterItn = params.Results.printitn;
stoptime = params.Results.stoptime;
stoptol = params.Results.stoptol;
% Extract the number of modes in tensor X.
N = ndims(X);
% If the initial guess has any rows of all zero elements, then modify
% so the row subproblem is not taking log(0). Values will be restored to
% zero later if the unfolded X for the row has no nonzeros.
for n = 1:N
rowsum = sum(Minit{n},2);
tmpIx = find(rowsum == 0);
if (isempty(tmpIx) == false)
Minit{n}(tmpIx,1) = 1.0e-8;
end
end
% Start with the initial guess, normalized using the vector L1 norm.
M = normalize(Minit,[],1);
% Sparse tensor flag affects how Pi and Phi are computed.
if isa(X,'sptensor')
isSparse = true;
else
isSparse = false;
end
% Initialize output arrays.
fnEvals = zeros(maxOuterIters,1);
kktViolations = -ones(maxOuterIters,1);
nInnerIters = zeros(maxOuterIters,1);
nzeros = zeros(maxOuterIters,1);
times = zeros(maxOuterIters,1);
if (printOuterItn > 0)
fprintf('\nCP_PDNR (alternating Poisson regression using damped Newton)\n');
end
dispLineWarn = (printInnerItn > 0);
% Start the wall clock timer.
tic;
if (isSparse && precomputeSparseIndices)
% Precompute sparse index sets for all the row subproblems.
% Takes more memory but can cut execution time significantly in some cases.
if (printOuterItn > 0)
fprintf(' Precomputing sparse index sets...');
end
sparseIx = cell(N);
for n = 1:N
num_rows = size(M{n},1);
sparseIx{n} = cell(num_rows,1);
for jj = 1:num_rows
sparseIx{n}{jj} = find(X.subs(:,n) == jj);
end
end
if (printOuterItn > 0)
fprintf('done\n');
end
end
e_vec = ones(1,R);
rowsubprobStopTol = stoptol;
%% Main Loop: Iterate until convergence or a max threshold is reached.
for iter = 1:maxOuterIters
isConverged = true;
kktModeViolations = zeros(N,1);
countInnerIters = zeros(1,N);
% Alternate thru each factor matrix, A_1, A_2, ... , A_N.
for n = 1:N
% Shift the weight from lambda to mode n.
M = redistribute(M,n);
% Calculate Khatri-Rhao product of all matrices but the n-th.
if (isSparse == false)
% Data is not a sparse tensor.
Pi = tt_calcpi_prowsubprob(X, isSparse, M, R, n, N, []);
X_mat = double(tenmat(X,n));
end
num_rows = size(M{n},1);
isRowNOTconverged = zeros(1,num_rows);
% Loop over the row subproblems in mode n.
for jj = 1:num_rows
% Initialize the damped Hessian parameter for the row subproblem.
mu = mu0;
% Get data values for row jj of matricized mode n.
if (isSparse)
% Data is a sparse tensor.
if (precomputeSparseIndices == false)
sparse_indices = find(X.subs(:,n) == jj);
else
sparse_indices = sparseIx{n}{jj};
end
if (isempty(sparse_indices))
% The row jj of matricized tensor X in mode n is empty.
M{n}(jj,:) = 0;
continue
end
x_row = X.vals(sparse_indices);
% Calculate just the columns of Pi needed for this row.
Pi = tt_calcpi_prowsubprob(X, isSparse, M, ...
R, n, N, sparse_indices);
else
x_row = X_mat(jj,:);
end
% Get current values of the row subproblem variables.
m_row = M{n}(jj,:);
% Iteratively solve the row subproblem with projected Newton steps.
innerIterMaximum = maxInnerIters;
if (inexactNewton && (iter == 1))
innerIterMaximum = 2;
end
for i = 1:innerIterMaximum
% Calculate the gradient.
[phi_row, ups_row] ...
= calc_partials(isSparse, Pi, epsDivZero, x_row, m_row);
gradM = (e_vec - phi_row)';
% Compute the row subproblem kkt_violation.
% Experiments in the original paper used this:
%kkt_violation = norm(abs(min(m_row,gradM')),2);
% Now we use | KKT |_inf:
kkt_violation = max(abs(min(m_row,gradM')));
% Report largest row subproblem initial violation.
if ((i == 1) && (kkt_violation > kktModeViolations(n)))
kktModeViolations(n) = kkt_violation;
end
if (mod(i, printInnerItn) == 0)
fprintf(' Mode = %1d, Row = %d, InnerIt = %d', ...
n, jj, i);
if (i == 1)
fprintf(', RowKKT = %.2e\n', kkt_violation);
else
fprintf(', RowKKT = %.2e, RowObj = %.4e\n', ...
kkt_violation, -f_new);
end
end
% Check for row subproblem convergence.
if (kkt_violation < rowsubprobStopTol)
break;
else
% Not converged, so m_row will be modified.
isRowNOTconverged(jj) = 1;
end
% Calculate the search direction.
[search_dir, predicted_red] ...
= getSearchDirPdnr(Pi, ups_row, R, gradM, m_row, mu, epsActSet);
% Perform a projected linesearch and update variables.
% Start from a unit step length, decrease by 1/2, stop with
% sufficient decrease of 1.0e-4 or at most 10 steps.
[m_rowNEW, f_old, f_unit, f_new, num_evals] ...
= tt_linesearch_prowsubprob(search_dir', gradM', m_row, ...
1, 1/2, 10, 1.0e-4, ...
isSparse, x_row, Pi, ...
phi_row, dispLineWarn);
fnEvals(iter) = fnEvals(iter) + num_evals;
m_row = m_rowNEW;
% Update damping parameter mu based on the unit step length,
% which is returned in f_unit.
actual_red = f_old - f_unit;
rho = actual_red / (-predicted_red);
if (predicted_red == 0)
mu = 10 * mu;
elseif (rho < 1/4)
mu = (7/2) * mu;
elseif (rho > 3/4)
mu = (2/7) * mu;
end
end
M{n}(jj,:) = m_row;
countInnerIters(n) = countInnerIters(n) + i;
end
% Test if all row subproblems have converged, which means that
% no variables in this mode were changed.
if (sum(isRowNOTconverged) ~= 0)
isConverged = false;
end
% Shift weight from mode n back to lambda.
M = normalize(M,[],1,n);
% Total number of inner iterations for a given outer iteration,
% totalled across all modes and all row subproblems in each mode.
nInnerIters(iter) = nInnerIters(iter) + countInnerIters(n);
end
% Save output items for the outer iteration.
num_zero = 0;
for n = 1:N
num_zero = num_zero + nnz(find(M{n} == 0.0));
end
nzeros(iter) = num_zero;
kktViolations(iter) = max(kktModeViolations);
if (inexactNewton)
rowsubprobStopTol = max(stoptol, kktViolations(iter) / 100.0);
end
% Print outer iteration status.
if (mod(iter,printOuterItn) == 0)
fprintf('%4d. Ttl Inner Its: %d, KKT viol = %.2e, obj = %.8e, nz: %d\n', ...
iter, nInnerIters(iter), kktViolations(iter), tt_loglikelihood(X,M), ...
num_zero);
end
times(iter) = toc;
% Check for convergence
if (isConverged && (inexactNewton == false))
break;
end
if (isConverged && (inexactNewton == true) && (rowsubprobStopTol <= stoptol))
break;
end
if (times(iter) > stoptime)
fprintf('Exiting because time limit exceeded\n');
break;
end
end
t_stop = toc;
%% Clean up final result and set output items.
M = normalize(M,'sort',1);
loglike = tt_loglikelihood(X,M);
if (printOuterItn > 0)
% For legacy reasons, compute "fit", the fraction explained by the model.
% Fit is in the range [0,1], with 1 being the best fit.
normX = norm(X);
normresidual = sqrt( normX^2 + norm(M)^2 - 2 * innerprod(X,M) );
fit = 1 - (normresidual / normX);
fprintf('===========================================\n');
fprintf(' Final log-likelihood = %e \n', loglike);
fprintf(' Final least squares fit = %e \n', fit);
fprintf(' Final KKT violation = %7.7e\n', kktViolations(iter));
fprintf(' Total inner iterations = %d\n', sum(nInnerIters));
fprintf(' Total execution time = %.2f secs\n', t_stop);
end
out = struct;
out.params = params.Results;
out.obj = loglike;
out.kktViolations = kktViolations(1:iter);