-
Notifications
You must be signed in to change notification settings - Fork 3
/
blurosy.m
591 lines (536 loc) · 22 KB
/
blurosy.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
function varargout=blurosy(th,params,xver,method,tsto,dth)
% [SbarordSbardth,k,tyy,CyyordCyydth]=blurosy(th,params,xver,method,tsto,dth)
%
% Wavenumber blurring of a univariate Matern spectral density (or its spectral
% derivative) with the periodogram of a spatial taper. The result is the expected
% periodogram (or a term that enters the variance calculation).
%
% Exact, fast, explicit way, no convolutional grid refinement (unlike BLUROS).
% Also unlike BLUROS, this IS a stand-alone code: it is not blurring an input
% parameterized spectral density but rather producing a blurred spectral density
% directly from input spectral parameters, through Fourier transformation of the
% product of the Matern correlation function with the autocorrelation of the
% taper window function.
%
% Equation numbers refer to Guillaumin et al., 2022, doi: 10.1111/rssb.12539
%
% INPUT:
%
% th The spectral parameter vector with elements:
% th(1)=s2 The first Matern parameter, aka sigma^2
% th(2)=nu The second Matern parameter
% th(3)=rho The third Matern parameter
% params Parameters of this experiment, the ones that are needed are:
% dydx sampling interval in the y and x directions [m m]
% NyNx number of samples in the y and x directions
% blurs -1 as appropriate for this procedure, any other errors
% taper 0 there is no taper near of far
% 1 it's a unit taper, implicitly
% OR an appropriately sized taper with explicit values
% (1 is yes and 0 is no and everything in between)
% xver 1 or 2 extra verification via BLURCHECK and alternative computations
% 0 No checking at all
% method 'ef' exact, efficient and fast [default]
% 'efs' exact, efficient and exploiting symmetry
% tsto An extra parameter slot to pass onto demo2
% dth 1, 2, or 3 specifies which element of th gets differentiated, to serve as
% an input to the variance calculation, which requires a "blurred" version of
% mAosl, if you will, but there is no sense doing it there.
%
% OUTPUT:
%
% SbarordSbardth The blurred spectral matrix or its derivative, with the unwrapped
% requested dimension as identified by the input 'params' (wrap with V2S)
% k The wavenumber matrix (the norm of the wave vectors), unwrapped
% tyy The autocorrelation of the spatial taper... which you
% may never need explicitly, used in SIMULOSL and LOGLIOSL
% CyyordCyydth The equivalent modified Matern correlation, may never need it explicitly
%
% SEE ALSO:
%
% SIMULOSL, BLUROS, MATERNOSP, BLURCHECK
%
% EXAMPLE:
%
% BLUROSY('demo1',pp,bb) compares against BLUROS where
% pp A square matrix size (one number)
% bb A MATERNOSP blurring densification (one number)
%
% BLUROSY('demo2',pp,nn,mm,tt) Boxcar/Ukraine/France example
% pp A matrix size (two numbers) [defaulted]
% nn A number of iterations to average over [defaulted]
% mm A method, either 'ef' or 'efs' [defaulted]
% tt 1 Boxcar 2 France 3 Ukraine
%
% BLUROSY('demo3') % should produce no output
%
% Last modified by arthur.guillaumin.14-at-ucl.ac.uk, 10/15/2017
% Last modified by olwalbert-at-princeton.edu, 06/17/2024
% Last modified by fjsimons-at-alum.mit.edu, 06/17/2024
if ~isstr(th)
if params.blurs>=0 & ~isinf(params.blurs)
error('Are you sure you should be running BLUROSY, not BLUROS?')
end
if isinf(params.blurs)
error('Are you sure you should be running BLUROSY, not MATERNOSY?')
end
% Defaults
defval('xver',1)
defval('method','ef')
defval('dth',[])
% Target dimensions, the original ones
NyNx=params.NyNx;
dydx=params.dydx;
switch method
case 'ef'
% Generates a double grid from which we subsample
% Fully exact and not particularly fast, still much faster than BLUROS
% Here are the full lags
ydim=[-NyNx(1):NyNx(1)-1]';
xdim=[-NyNx(2):NyNx(2)-1] ;
% Here is the Matern spatial covariance on the double distance grid,
% multiplied by the spatial taper in a way that its Fourier
% transform can be the convolution of the spectral density with the
% spectral density of the taper, i.e. the expected periodogram
[Cyy,tyy]=spatmat(ydim,xdim,th,params,xver,dth);
% http://blogs.mathworks.com/steve/2010/07/16/complex-surprises-from-fft/
% Here is the blurred covariance on the 'double' grid
Hh=fftshift(realize(fft2(ifftshift(Cyy))));
% Play with a culled DFTMTX? Rather now subsample to the 'complete' grid
Hh=Hh(1+mod(NyNx(1),2):2:end,1+mod(NyNx(2),2):2:end);
case 'efs'
% Generates a sample-size grid by working from a quarter, rest symmetric
% Fully exact and trying to be faster for advanced symmetry in the covariance
% Here are the partial lags
ydim=[0:NyNx(1)-1]';
xdim=[0:NyNx(2)-1] ;
% Here is the Matern spatial covariance on the quarter distance grid, see above
[Cyy,tyy]=spatmat(ydim,xdim,th,params,xver,dth);
% Exploit the symmetry just a tad, which allows us to work with smaller matrices
q1=fft2(Cyy);
q4=q1+[q1(:,1) fliplr(q1(:,2:end))];
% Here is the blurred covariance on the 'complete' grid
Hh=fftshift(2*real(q4-repmat(fft(Cyy(:,1)),1,NyNx(2)))...
-repmat(2*real(fft(Cyy(1,1:end))),NyNx(1),1)...
+Cyy(1,1));
if nargout>2
% If you ever wanted tyy/Cyy to come out you'll need to unquarter it
tyy=[fliplr(tyy(:,2:end)) tyy]; tyy=[flipud(tyy(2:end,:)) ; tyy];
tyy=[zeros(size(tyy,1)+1,1) [zeros(1,size(tyy,2)) ; tyy]];
Cyy=[fliplr(Cyy(:,2:end)) Cyy]; Cyy=[flipud(Cyy(2:end,:)) ; Cyy];
Cyy=[zeros(size(Cyy,1)+1,1) [zeros(1,size(Cyy,2)) ; Cyy]];
end
end
% Normalize and vectorize
Sbar=Hh(:)*prod(dydx)/(2*pi)^2;
% Should check positivity always!
if any(Sbar<0); keyboard; end
% Check Hermiticity of the results
if xver==1 || xver==2
blurcheck(Sbar,params)
if strcmp('method','ef')
hermcheck(tyy)
hermcheck(Cyy)
end
end
% Produce the unwrapped wavenumbers if you've requested them to be output
if nargout>1
k=knums(params);
k=k(:);
else
k=[];
end
% Optional output
varns={Sbar,k,tyy,Cyy};
varargout=varns(1:nargout);
elseif strcmp(th,'demo1')
% This case is handled by BLUROS
% Remember the prompt parameter names are no longer what they seem, the first
% one is matrix size and the second one blurring refinement density den
bluros('demo1',params,xver)
elseif strcmp(th,'demo2')
% Now here we will test that the blurred spectrogram is the average
% periodogram of the periodogram of the sample generated by the spectral
% density. We welcome variability around it, as this is in line with
% the uncorrelated-wavenumber approximation of the debiased Whittle
% approach, but we do fear systematic offsets.
% Number of iterations
defval('xver',100);
% Method of computation
defval('method','efs')
% Type of taper, e.g. boxcar, France, Ukraine
defval('tsto',1)
% Field size
defval('params',[188 233]+randi(20,[1 2]))
% BEGIN Figure for boxcar paper
params=[33 33]*3;
% Some Matern paramters
th=[1 2.5 1e3];
% Some number of iterations
xver=100;
% Taper type
tsto=1;
% END Figure for boxcar paper
% Some combinations - SIMULOSL
p.NyNx=[params(1) params(2)];
p.dydx=1e3*[1 1];
% Compare the average periodogram with the blurred spectral density
% What kind of a test are we running? Boxcar, or France/Ukraine?
switch tsto
case 1
% Boxcar is 0 or 1 or ones
p.taper=1;
p.mask='boxcar';
% The mask parameter is irrelevant, though it might be the anti-taper
case 2
% Here's another one
p.mask='france';
case 3
% Here's another one
p.mask='ukraine';
end
switch tsto
case {2,3}
% Generate "mask" only, never mind what the data will be
[~,~,I]=maskit(rand(p.NyNx),p);
% Keep the mask as a taper or the taper as mask, for illustration only
p.taper=I;
end
% Calculate expected periodogram, i.e. the appropriately blurred likelihood
% Use the exact method via BLUROSY, force p.blurs=-1
p.blurs=-1;
% Just do the real xver=1 explicitly here
Sbar=blurosy(th,p,1,method);
% Then for what comes next, to simulate data using SGP, force p.blurs=Inf
p.blurs=Inf;
% One random one from the sequence will be shown
randix=randi(xver);
% We be collecting the average periodogram to compare with the expected periodogram
Sbb=0;
% The third input in the demo was the number of iterations, fake-called 'xver'
for index=1:xver
% Simulate sequentially, collect the expected periodogram, and
% make the average periodogram
[Hx,th0,p,k,Hk,Sb,Lb,gane,miy]=simulosl(th,p);
% Are none of them really bad?
if max(Sb./Sbar)>20; warning('ratio seems bad') ; end
if index==randix;
% Keep the special ones for plotting later
Hxx=Hx; Sbx=Sb;
% This is the special one that we shall plot later on
Xk=Sb./Sbar; varibal='X';
end
% Collect the average, watch the growth later
Sbb=Sbb+Sb;
% Mean of the evolving ratio over the wavenumbers
m(index)=mean(Sbb/index./Sbar);
% Mean of the evolving standard deviation of the mean ratio over the wavenumbers
s(index)=std(Sbb/index./Sbar);
end
% Mean ratio over the realizations
Sbb=Sbb/xver;
% Number of standard deviations for the axis limits and color bars
nsig=2;
% First figure %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
figure(1)
clf
[ah,ha,H]=krijetem(subnum(2,3));
axes(ah(1))
imagefnan([1 1],p.NyNx([2 1]),v2s(Hxx,p),[],th(1)*[-3 3]); axis image ij
t(1)=title(sprintf('field # %i',randix));
axes(ah(2))
imagesc(log10(v2s(Sbx,p))); axis image
t(2)=title(sprintf('periodogram # %i',randix));
% The spatial domain taper
axes(ah(4))
if length(p.taper)==1; p.taper=ones(p.NyNx); end
imagesc(p.taper); axis image
t(4)=title('taper');
% The expectation of the periodogram
axes(ah(3))
imagesc(log10(v2s(Sbar,p))); axis image
t(3)=title(sprintf('expectation | %s [%g %g %gx]',...
'\theta =',th./[1 1 sqrt(prod(p.dydx))]));
% Then compare with the thing coming out of BLUROSY
axes(ah(6))
imagesc(log10(v2s(Sbb,p))); axis image
t(6)=title(sprintf('average | %i realizations',xver));
axes(ah(5))
imagesc(v2s(Sbb./Sbar,p)); axis image
t(5)=title(sprintf('(aver / expec), m %4.2f, s %4.2f',...
m(end),s(end)));
try
set(ah(5),'clim',m(end)+[-1 1]*nsig*s(end))
end
% Clean that sh*t up
[k,dci,dcn,kx,ky]=knums(p);
set(ah([1 4]),'xtick',unique([1 dci(2) p.NyNx(2)]),...
'ytick',unique([1 dci(1) p.NyNx(1)]));
set(ah([2 3 5 6]),'xtick',unique([1 dci(2) p.NyNx(2)]),...
'ytick',unique([1 dci(1) p.NyNx(1)]),...
'xticklabel',[-1 0 1],...
'yticklabel',[-1 0 1]);
longticks(ah)
set(ah([3 6]),'YAxisLocation','right')
% serre(H,1,'across')
serre(H',0.5,'down')
movev(t,-p.NyNx(1)/20)
figdisp([],sprintf('demo_2_%3.3i_%3.3i_%3.3i-%3.3i_fld',p.NyNx,tsto,xver),[],2)
% Second figure %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
figure(2)
clf
clf
[ah,ha,H]=krijetem(subnum(2,3));
df=2;
% Craft some labels
xll=[0 3*2*df];
xlls=[xll(1):df:xll(2)];
xstr=sprintf('quadratic residual %s',varibal);
xstr2=sprintf('quadratic residual 2%s',varibal);
cax=[0 3*df];
ly=p.NyNx(1);
lx=p.NyNx(2);
ylis=unique([1 dci(2) p.NyNx(2)]);
xlis=unique([1 dci(1) p.NyNx(1)]);
cborien='vert';
axes(ah(1))
% Literally follow MLECHIPLOS, to a point, spin off function later on
[bdens,c]=hist(2*Xk,5*round(log(prod(p.NyNx))+1));
% Plot the histogram as a bar graph
bdens=bdens/indeks(diff(c),1)/prod(p.NyNx);
bb=bar(c,bdens,1);
% Each of the below should be df/2
t(1)=title(sprintf('m(%s) = %5.3f v(%s) = %5.3f',...
varibal,nanmean(Xk),...
varibal,nanvar(Xk)));
set(bb,'FaceC',grey)
hold on
% Plot the ideal chi-squared distribution
refs=linspace(0,max(2*Xk),100);
hold on
plot(refs,chi2pdf(refs,df),'Linew',0.5,'Color','k')
hold off
longticks(ah(1))
maxi=0.25;
ylls=[0:0.1:maxi*(4/df)];
ylim([0 maxi*(4/df)])
xlim(xll)
set(ah(1),'xtick',xlls)
xl(1)=xlabel(xstr2);
yl(1)=ylabel('probability density');
axis square
movev(t(1),maxi*(4/df)/20)
axes(ah(2))
% This should be a straight line folks since the ratio is 1/2 chi^2_2
h=qqplot(2*Xk,makedist('gamma','a',df/2,'b',2)); axis equal;
axis image; box on
longticks(ah(2))
set(h(1),'MarkerEdge','k')
set(h(3),'LineStyle','-','Color',grey)
% Plot the one-to-one line
hold on
xh=get(h,'Xdata'); xh=xh{1};
yh=get(h,'ydata'); yh=yh{1};
qq0=plot([0 xll(2)],[0 xll(2)],'k');
delete(h(1))
% Replot for more control
plot(xh,yh,'LineStyle','none','Marker','+','MarkerFaceColor','k',...
'MarkerEdgeColor','k','MarkerSize',2);
top(h(3),ah(2))
delete(get(ah(2),'ylabel'));
delete(get(ah(2),'title'));
xlim(xll); ylim(xll)
xl(2)=xlabel(sprintf('predicted 2%s',varibal));
yl(2)=ylabel(sprintf('observed 2%s',varibal));
set(gca,'xtick',xlls,'ytick',xlls)
allg=~isinf(Xk); sum(allg);
% Test for departure of chi-squaredness via the "magic" parameter which
% This is the same as what comes out of LOGLIOSL etc
magx=nanmean([Xk-df/2].^2);
neem='\xi'; neem='s_X^2';
% Do the test whether you accept this as a good fit, or not
vr=8/length(k(~~k));
[a,b,c]=normtest(magx,1,vr,0.05);
% disp(sprintf('NORMTEST %i %5.3f %i',a,b,round(c)))
if a==0; stp='accept'; else stp='reject'; end
t(2)=title(sprintf('%s = %5.3f 8/K = %5.3f %s p = %5.2f',...
neem,magx,vr,stp,b));
movev(t(2),xll(2)/40)
axes(ah(3))
imagefnan([xlis(1) ylis(end)],[xlis(end) ylis(1)],...
v2s(Xk,p),'gray',cax,[],1); axis image square
set(ah(3),'xtick',xlis,'XtickLabel',[-1 0 1],...
'ytick',ylis,'YtickLabel',[-1 0 1])
longticks(ah(3))
xl(3)=xlabel('relative wavenumber');
yl(3)=ylabel('relative wavenumber');
t(3)=title(sprintf('%s | %i x %i | %s = [%g %g %gx]',p.mask,p.NyNx,'\theta',...
th./[1 1 sqrt(prod(p.dydx))]));
movev(t(3),ylis(end)/20)
delete(ah([4:6]))
% Cosmetics from EGGERS6
axes(ah(3))
[cb,xcb]=addcb(cborien,cax,cax,'gray',df,1);
axes(cb)
set(xcb,'string',xstr)
set(cb,'YAxisLocation','right')
set(cb,'position',...
[getpos(ah(3),1)+getpos(ah(3),3)*1.1 getpos(ah(3),2) getpos(cb,3) getpos(ah(3),4)])
shrink(cb,1,1.15)
% Control all axes and font sizes at the same time
set(ah(1:3),'FontSize',8)
set(t,'FontSize',8-1)
figdisp([],sprintf('demo_2_%3.3i_%3.3i_%3.3i-%3.3i_chi',p.NyNx,tsto,xver),[],2)
% Third figure %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
figure(3)
clf
H=errorbar(1:xver,m,-nsig*s,nsig*s);
try
set(getkids(H,1),'LineWidth',0.5)
set(getkids(H,2),'Color',grey);
catch
H.LineWidth=0.5;
H.Color=grey;
end
longticks(gca,2)
axis tight; grid on; xlabel('sample size');
ylabel(sprintf('average / expected periodogram | %i sigma',nsig))
xlim([0 xver+1])
set(gca,'xtick',[1 10:10:xver])
shrink(gca,1,1.5)
ylim([min(m-s*2) max(m+2*s)]*1.1)
hold on
plot(1:xver,m,'LineWidth',1,'Color','b')
hold off
tt=title(sprintf('%s | %i x %i | %s = [%g %g %gx]',p.mask,p.NyNx,'\theta',...
th./[1 1 sqrt(prod(p.dydx))]));
movev(tt,range(ylim)/30)
figdisp([],sprintf('demo_2_%3.3i_%3.3i_%3.3i-%3.3i_std',p.NyNx,tsto,xver),[],2)
elseif strcmp(th,'demo3')
% Simulate some random data with default p.blurs=Inf and p.taper=0
[H,th,p]=simulosl;
% Reset to implicit treatment of unit taper
p.blurs=-1; p.taper=1;
% This function produces the blurred spectral densities, different methods
[Sbar1,k1,tyy1,Cyy1]=blurosy(th,p,1,'ef');
[Sbar2,k2,tyy2,Cyy2]=blurosy(th,p,1,'efs');
% Reset to explicit treatment of unit taper
p.taper=ones(p.NyNx);
[Sbar3,k3,tyy3,Cyy3]=blurosy(th,p,1,'ef');
[Sbar4,k4,tyy4,Cyy4]=blurosy(th,p,1,'efs');
% All of these should be virtually identical
% Note the Sbar are large so the tolerance is high
tolex=-(log10(max(Sbar1))-11);
diferm(Sbar1,Sbar2,tolex);
diferm(tyy1,tyy2); diferm(Cyy1,Cyy2);
diferm(Sbar1,Sbar3,tolex);
diferm(tyy1,tyy3); diferm(Cyy1,Cyy3);
diferm(Sbar1,Sbar4,tolex);
diferm(tyy1,tyy4); diferm(Cyy1,Cyy4);
diferm(Sbar2,Sbar3,tolex);
diferm(tyy2,tyy4); diferm(Cyy2,Cyy3);
diferm(Sbar2,Sbar4,tolex);
diferm(tyy2,tyy4); diferm(Cyy2,Cyy4);
diferm(Sbar3,Sbar4,tolex);
diferm(tyy3,tyy4); diferm(Cyy3,Cyy4);
% subplot(221); imagesc(log10(v2s(Sbar1,p))); axis square
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function [Cyy,t]=spatmat(ydim,xdim,th,params,xver,dth)
% [Cyy,t]=spatmat(ydim,xdim,th,params,xver,dth)
%
% Returns the modified spatial covariance whose Fourier transform is the blurred
% spectrum after spatial data tapering, i.e. the expected periodogram. No use to
% return the autocovariance of the applied spatial taper which is an essential
% part of this operation... never need it explicitly, used in SIMULOSL and
% LOGLIOSL via the intermediary of BLUROSY.
% Dimensions of the original grid
NyNx=params.NyNx;
dydx=params.dydx;
if ~isfield(params,'taper')
params.taper=[];
end
% Specify the spatial taper EXPLICITLY
if prod(size(params.taper))>1
% Now compare with the other mechanisms
% Completely general windows where 1 means you are taking a sample
% See Arthur's note for more general windows, use IFF2/FFT2 you need, see
% ~/POSTDOCS/ArthurGuillaumin/CodeArthur/NonParametricEstimation/Periodogram.m
% $MFILES/retired/QR?.m/map_*.m/whittle_*/expected_* etc
Tx=double(params.taper);
% If you are here with 'efs' the taper is explicit AND not
% symmetric, so must do something else
if all([length(ydim) length(xdim)]==NyNx)
% Produce the autocorrelation sequence eq. (12)
t=zeros(size(Tx));
% It's quite vital that these be colon ranges (faster) or (like
% here) ROW index vectors... mixing rows/columns won't work
for i=ydim(:)'+1
for j=xdim(:)'+1
% Vectorize? Check out XCORR2, that's also good
t(i,j)=sum(sum(Tx(1:NyNx(1)-i+1,1:NyNx(2)-j+1).*(conj(Tx(i:end,j:end)))));
end
end
else
% This also obviates the extra test below, really this should be the
% top way, but leave the explicit way for illustration - this is
% the long step
t=xcorr2(Tx);
% Add a row of zeros here
t=[zeros(size(t,1)+1,1) [zeros(1,size(t,2)) ; t]];
end
% Now normalize the cross-correlations at the end
t=t/sum(sum(Tx.^2));
% disp(sprintf('\nStill can optimize\n'))
% Here too should use FFT where we can, see COMPUTE_KERNELS
% internally and below, I would imagine that's just the same thing
% inside Arthur's PERIODOGRAM object
% kernel1 = ifft2(abs(fft2(Tx,2*p.NyNx(1),2*NyNx(2))).^2/prod(p.NyNx));
% kernel1 = kernel1(1:p.NyNx(1),1:p.NyNx(2));
% kernel2 = ifft2(abs(fft2(fliplr(Tx),2*p.NyNx(1),2*p.NyNx(2))).^2/prod(p.NyNx));
% kernel2 = kernel2(1:p.NyNx(1),1:p.NyNx(2));
% Then feed that into COMPUTE_EXPECTATION_FROM_KERNELS
if (xver==1 || xver==2) && all(size(t)==NyNx)
t3=xcorr2(Tx); t3=t3/sum(sum(Tx.^2));
if all(size(t)==NyNx)
% Then the doubling inside this block needs to be undone
t3=t3(NyNx(1):end,NyNx(2):end);
end
diferm(t,t3);
end
else
% Specify the spatial taper IMPLICITLY, taper is just a single number, could
% be 0 or 1 both telling us "not" spatially tapered which is of course in
% essence the same as a "unit" spatial taper taper operation. The triangle
% functions are the normalized autocorrelations of the unit window functions,
% i.e. c_{g,n}(u) of (12)-(13) in Guillaumin (2022), doi: 10.1111/rssb.12539
triy=1-abs(ydim)/NyNx(1);
trix=1-abs(xdim)/NyNx(2);
% Here is the gridded triangle for this case
t=bsxfun(@times,triy,trix);
if xver==1 || xver==2
% Do form the taper explicitly after all, normalize ahead of time
Tx=ones(NyNx)/sqrt(prod(NyNx));
% Need to cut one off Arthur says, possibly need to
% re-re-visit these even/odd comparisons in BLUROS, if it ever
% gets to that point; currently the comparison is favorable
t2=fftshift(ifft2(abs(fft2(Tx,2*size(Tx,1)-1,2*size(Tx,2)-1)).^2));
% Fix the rim by adding zeroes top and left
t2=[zeros(size(t2,1)+1,1) [zeros(1,size(t2,2)) ; t2]];
% Check the difference between these two implementations,
% all checked for even/odd/method combinations on 2/24/2023
if all(size(t)==NyNx)
% Then the doubling inside this block needs to be undone
t2=t2(NyNx(1)+1:end,NyNx(2)+1:end);
end
diferm(t,t2,9);
end
end
% Here is the distance grid, whose size depends on the input
y=sqrt(bsxfun(@plus,[ydim*dydx(1)].^2,[xdim*dydx(2)].^2));
% The modified spatial covariance
Cyy=maternosy(y,th,[],dth).*t;
% Remind me: Arthur: diag(U^T * Cyy * U) where U is the DFMTX is the
% expected periodogram, the diagonal variance, but to get the variance of
% the gradient we need the whole thing, the covariance, (U^T * Cyy * U) of
% course in two dimensions, so properly arranged.