-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathsvdslep1.m
282 lines (243 loc) · 7.49 KB
/
svdslep1.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
function varargout=svdslep1(N,NW,K,method,imp)
% [E,V,SE,ngro,]=SVDSLEP1(N,NW,K,method,imp)
%
% Explicit diagonalization of the Slepian concentration operators
% in one linear dimension.
%
% INPUT:
%
% N The length of the sequence [default: 2^8]
% NW Half the Shannon number (time-bandwidth product) [default: 10]
% K The requested number of eigenvectors [default: 2NW-1]
% method 1 Using EIGS on the Hermitian form [default]
% 2 Using EIG on the Hermitian form
% 3 Using SVDS on the projection operator
% 4 Using SVD on the projection operator
% 5 Using the power method for the largest eigenvalue
% imp 1 Uses the implicit operator method [default]
% 0 Uses the sparse explicit operator method [slow]
%
% OUTPUT:
%
% E The eigenfunctions of the concentration problem
% V The eigenvalues of the concentration problem
% SE The power spectrum of the eigenfunctions
% ngro The "growth" factor used
%
% EXAMPLE:
%
% svdslep1('demo1')
%
% Last modified by fjsimons-at-alum.mit.edu, 02/19/2020
% See some notes in RB VIII p 90
% Default values
defval('N',2^10)
if ~isstr(N)
defval('NW',10)
defval('K',2*NW-1)
defval('method',1)
defval('ngro',256);
defval('imp',1);
% This one for "excessive" verification
defval('xver',0)
% Make a larger domain
Nd=N*ngro;
if ngro>32 || N>256
% Force to the implicit method or your machine will die
imp=1;
disp('Method reset to IMPLICIT')
end
if imp==1
disp(sprintf('Implicit method embedded %ix\n',ngro))
nonz=(Nd-N)/2+[1:N];
P =@(x) proj(x,nonz);
Q =@(x) fft(x);
Qi=@(x) ifft(x);
nonz=[1:ngro*NW+1 Nd-ngro*NW+1:Nd];
L =@(y) proj(y,nonz);
H =@(x) P(Qi(L(Q(P(x)))));
% Acknowledge that H is complex (though it is symmetric)
OPTS.isreal=false;
OPTS.disp=0;
% Remember to specify the output size
[E,V]=eigs(H,Nd,K,'LR',OPTS);
[V,i]=sort(diag(V),'descend');
E=E(:,i); V=V(1:K); E=E(:,1:K);
elseif imp==0
disp(sprintf('Explicit, sparse method embedded %ix\n',ngro))
% The spatial projection operator
P=diag([zeros(1,(Nd-N)/2) ones(1,N) zeros(1,(Nd-N)/2)]);
% The Fourier transform operator
Q=dftmtx(Nd)/sqrt(Nd);
% The bandlimiting operator
L=diag([ones(1,ngro*NW+1) zeros(1,Nd-2*ngro*NW-1) ones(1,ngro*NW)]);
% Make them all sparse for speed - could start from it altogether
P=sparse(P);
Q=sparse(Q);
L=sparse(L);
if xver==1
% Check that Q is a unitary matrix
difer(Q'*Q-diag(diag(Q'*Q)))
difer(Q'*Q-eye(size(Q)))
% Check that H is a Hermitian matrix
difer(H'-H)
% Check that H is a positive definite matrix
[R,p]=chol(H); isposdef=p==0; difer(isposdef-1)
% Check that Q does what I think it does
f=rand(Nd,1);
% See the normalization of FFT (1) and IFFT (1/N)
difer(fft(f)/sqrt(Nd)-Q*f);
end
% The operator whose singular functions we want
A=L*Q*P;
% The operator whose eigenfunctions we want
H=P*Q'*A;
% The singular values of A are the eigenvalues of AA'. The singular
% values of the concentration/projection operator are the eigenvalues of
% the "squared" projection operator.
% The eigenvector decomposition
switch method
case 1
disp('Using EIGS')
OPTS.disp=0;
[E,V]=eigs(H,K,'LR',OPTS);
[V,i]=sort(diag(V),'descend');
E=E(:,i); V=V(1:K); E=E(:,1:K);
case 2
disp('Using EIG')
[E,V]=eig(H,'nobalance');
[V,i]=sort(diag(V),'descend');
E=E(:,i); V=V(1:K); E=E(:,1:K);
case 3
disp('Using SVDS')
[U,V,E]=svds(A,K);
[V,i]=sort(diag(V).^2,'descend');
E=E(:,i); V=V(1:K); E=E(:,1:K);
case 4
disp('Using SVD')
[U,V,E]=svd(A);
[V,i]=sort(diag(V).^2,'descend');
E=E(:,i); V=V(1:K); E=E(:,1:K);
case 5
disp('Using the power method for the first one only')
Vi=1; Vj=2;
E=rand(ngro*N,1);
E=E/sqrt(E'*E);
while abs(Vj-Vi)>1e-6
Vi=Vj;
Vj=E;
E=H*E;
Vj=real(E'*Vj);
E=real(E/sqrt(E'*E));
plot(E);
pause(0.1)
end
V=Vj;
end
end
% Define some kind of tolerance level
tol=100*eps;
% Make them real as we know they should be
if any(imag(V)>tol)
error('Complex eigenvalues');
else
V=real(V);
% Check imaginary part of the "good" eigenfunctions
disp(sprintf('mean(abs(imag(E))) = %8.3e out of %8.3e\n',...
mean(mean(abs(imag(E(:,V>tol))))),...
mean(mean(abs(E(:,V>tol))))))
% Take out only the central part
E=E((Nd-N)/2+1:(Nd-N)/2+N,:);
% Note that they were normalized in the complex plane
E=real(E); E=E./repmat(diag(sqrt(E'*E))',size(E,1),1);
% Cannot do abs(E).*sign(E) because the sign is bad for complex
% But either real or imag seems to work... they more or less scale!
end
if xver==1
% Check the orthogonality for the "good" eigenfunctions
ortho=E(:,V>tol)'*E(:,V>tol);
difer(ortho-diag(diag(ortho)))
difer(diag(ortho)-1)
end
if nargout>2
% Get the power spectrum
SE=abs(fft(E,8*length(E))).^2;
else
SE=NaN;
end
% Fix the first bit to be an upswing
for index=1:size(E,2)
E(:,index)=E(:,index)*sign(E(2,index));
end
% Output
varns={E,V,SE,ngro};
varargout=varns(1:nargout);
elseif strcmp(N,'demo1')
% So you can say svdslep1('demo1',1)
defval('NW',[])
imp=NW;
defval('imp',1)
N=2^(3+ceil(rand*5)); NW=ceil(rand*5);
N=2^8; NW=3;
method=1;
disp(sprintf('\nN = %i , NW = %i\n',N,NW))
[E,V,SE,ngro]=svdslep1(N,NW,max(3,2*NW),method,imp);
% Compare with Matlab's built-in functions
[E2,V2]=dpss(N,NW,max(3,2*NW));
clf
[ah,ha]=krijetem(subnum(2,3));
% Make the frequency axis for the positive frequencies only
[fax,selekt]=fftaxis1D(E(:,1),length(SE),N-1);
SE=SE(selekt,:);
% Another, identical, more recently revised, way
% [K,kx]=knum2([2 N],[2 N]-1);
% fx=-fliplr(indeks(kx,selekt)/2/pi);
yls=[-100 2];
for ind=1:3
SE2=indeks(abs(fft(E2(:,ind),8*length(E2)).^2),selekt);
axes(ah(ind))
plot(E(:,ind),'b','LineWidth',2); hold on
plot(E2(:,ind),'y')
disp(sprintf('Eigenfunction %i: rms %8.3f%s',ind,...
rms(E(:,ind)-E2(:,ind))/rms(E(:,ind))*100,'%'))
yl(ind)=ylabel('Slepian eigenfunction');
set(ah(ind),'xtick',[1 N],'xlim',[1 N])
xlabel(sprintf('%s = %12.9f (EIGS)','\lambda',V(ind)))
axes(ha(2*ind))
plot(N*fax,decibel(SE(:,ind)),'b','LineWidth',2); hold on
plot(N*fax,decibel(SE2),'y');
plot(N*[NW/N NW/N],yls,'k-');
plot(-N*[NW/N NW/N],yls,'k-');
hold off
%xl(3+ind)=xlabel('frequency (sample)^{-1}');
xl(3+ind)=xlabel('frequency x sample size');
yl(3+ind)=ylabel('power spectral density (dB)');
title(sprintf('%s = %12.9f (DPSS)','\lambda',V2(ind)))
% Show the entire range
too=8*1.95;
exes=[0 1/too];
exes=[-0.1 20];
set(ha(2*ind),'xtick',round(exes*10)/10,'xlim',minmax(exes),...
'xgrid','off','ylim',yls)
end
% Cosmetics
delete(yl([2 3 5 6]))
longticks(ah)
tt=supertit(ah(1:3),...
sprintf('SVDSLEP1 vs DPSS on a %ix domain',ngro),12);
movev(tt,0.15)
fig2print(gcf,'landscape')
% Is the ratio of the two things symmetric?
% for ind=1:size(E,2)
% a=E(:,ind)./E2(:,ind);
% rmsym(ind)=rms([a(1:length(a)/2)-flipud(a(length(a)/2+1:end))])./...
% rms(a)*100;
% end
% rmsym
figdisp([],'demo1')
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function Pv=proj(v,p)
% Projects the vector v on the indices p
Pv=zeros(size(v));
Pv(p)=v(p);