-
Notifications
You must be signed in to change notification settings - Fork 28
/
phaseStats.m
executable file
·1757 lines (1523 loc) · 78.7 KB
/
phaseStats.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
classdef phaseStats
% Phase statistics static class
methods (Static)
function out = variance(atm)
%% VARIANCE Phase variance
%
% out = phaseStats.variance(atm) computes the phase variance
% from an atmosphere object
%
% See also atmosphere
L0r0ratio= (atm.L0./atm.r0).^(5./3);
out = (24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6).*gamma(5./6)./(2.*pi.^(8./3))).*L0r0ratio;
layers = atm.layer;
out = sum([layers.fractionnalR0]).*out;
end
function out = covariance(rho,atm)
%% COVARIANCE Phase covariance
%
% out = phaseStats.covariance(rho,atm) computes the phase covariance from
% the baseline rho and an atmosphere object
%
% See also atmosphere
L0r0ratio= (atm.L0./atm.r0).^(5./3);
cst = (24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6)./(2.^(5./6).*pi.^(8./3))).*...
L0r0ratio;
out = ones(size(rho)).*(24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6).*gamma(5./6)./(2.*pi.^(8./3))).*L0r0ratio;
index = rho~=0;
u = 2.*pi.*rho(index)./atm.L0;
out(index) = cst.*u.^(5./6).*besselk(5./6,u);
layers = atm.layer;
out = sum([layers.fractionnalR0]).*out;
end
function out = covarianceXDerivative(rho,argRho,atm)
%% COVARIANCE Phase covariance
%
% out = phaseStats.covariance(rho,atm) computes the phase covariance from
% the baseline rho and an atmosphere object
%
% See also atmosphere
L0r0ratio= (atm.L0./atm.r0).^(5./3);
cst = (24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6)./(2.^(5./6).*pi.^(8./3))).*...
L0r0ratio;
out = zeros(size(rho));
index = rho~=0;
u = 2.*pi.*rho(index)./atm.L0;
out(index) = -cst.*(1/6).*u.^(-1./6).*(2*pi/atm.L0).*cos(argRho(index)).*...
(-5*besselk(5./6,u) + 3.*u.*(besselk(-1./6,u)+besselk(11./6,u)));
layers = atm.layer;
out = sum([layers.fractionnalR0]).*out;
end
function out = angularCovariance(theta,atm)
%% ANGULARCOVARIANCE Phase angular covariance
%
% out = phaseStats.angularCovariance(theta,atm) computes the
% phase angular covariance from the zenith angle theta and an
% atmosphere object
%
% See also atmosphere
out = zeros(size(theta));
for kLayer = 1:atm.nLayer
atmSlab = slab(atm,kLayer);
out = out + phaseStats.covariance(atmSlab.layer.altitude*tan(theta),atmSlab);
end
end
function out = angularStructureFunction(theta,atm)
%% ANGULARSTRUCTUREFUNCTION Phase angular structure function
%
% out = phaseStats.angularStructureFunction(rho,atm) computes the
% phase angular structure function from the zenith angle theta
% and an atmosphere object
%
% See also atmosphere
out = zeros(size(theta));
for kLayer = 1:atm.nLayer
atmSlab = slab(atm,kLayer);
out = out + 2.*( phaseStats.variance(atmSlab) - ...
phaseStats.covariance(atmSlab.layer.altitude*tan(theta),atmSlab) );
end
end
function out = temporalCovariance(tau,atm)
%% TEMPORALCOVARIANCE Phase temporal covariance
%
% out = phaseStats.temporalCovariance(rho,atm) computes the
% phase temporal covariance from the delay tau and an
% atmosphere object
%
% See also atmosphere
out = zeros(size(tau));
for kLayer = 1:atm.nLayer
atmSlab = slab(atm,kLayer);
out = out + phaseStats.covariance(atmSlab.layer.windSpeed*tau,atmSlab);
end
end
function out = temporalStructureFunction(tau,atm)
%% TEMPORALSTRUCTUREFUNCTION Phase temporal structure function
%
% out = phaseStats.temporalStructureFunction(rho,atm) computes
% the phase temporal structure function from the delay tau and
% an atmosphere object
%
% See also atmosphere
out = zeros(size(tau));
for kLayer = 1:atm.nLayer
atmSlab = slab(atm,kLayer);
out = out + 2.*( phaseStats.variance(atmSlab) - ...
phaseStats.covariance(atmSlab.layer.windSpeed*tau,atmSlab) );
end
end
function L2 = sparseInverseCovarianceMatrix(layerGrid,pitch,atm)
%% SPARSEINVERSECOVARIANCEMATRIX
%
% out = sparseInverseCovarianceMatrix(gridMask,pitch,atm)
% computes a sparse approximation of the inverse of the
% covariance matrix from the pupil grid mask and the atmosphere
%
% See also atmosphere
% nGrid = length(layerGrid);
% rho = [0 1 nGrid]*pitch;
% cov = phaseStats.covariance(atm,rho);
% cov = cov.*[-4 2 2];
% sigma = sum(cov);
% L2 = cell(1,nGrid);
% for kGrid=1:nGrid
% gridMask = layerGrid{kGrid};
% [n,m] = size(gridMask);
% if any([n,m])==1
% n = sqrt(numel(gridMask));
% end
% L = gallery('poisson',n);
% L(~gridMask,:) = [];
% L(:,~gridMask) = [];
% L = (L'*L);
% L2{kGrid} = L./(phaseStats.variance(slab(atm,kGrid))*sum(L(:)));
% end
nGrid = length(layerGrid);
L2 = cell(1,nGrid);
for kGrid=1:nGrid
gridMask = layerGrid{kGrid};
n = sqrt(numel(gridMask));
e = ones(n^2,1);
% e(~gridMask) = [];
L = spdiags([e e -4*e e e],[-n -1 0 1 n],n^2,n^2);
L(~gridMask,:) = [];
L(:,~gridMask) = [];
L2{kGrid} = L'*L;
[i,j] = find(L2{kGrid});
[x,y] = meshgrid( (0:n-1)*pitch );
z = complex(x(:),y(:));
psCov = phaseStats.covariance(abs(z(i)-z(j)),slab(atm,kGrid));
% psCov = psCov*atm.la
alpha2 = sum(nonzeros(L2{kGrid}).*psCov);
L2{kGrid} = L2{kGrid}*n^2/alpha2;
end
end
function out = structureFunction(rho,atm)
%% STRUCTUREFUNCTION Phase structure function
%
% out = phaseStats.structureFunction(rho,atm) computes the
% phase structure function from the baseline rho and an
% atmosphere object
%
% See also atmosphere
if isinf(atm.L0)
out = zeros(size(rho));
index = rho~=0;
out(index) = 2.*(24.*gamma(6./5)./5).^(5./6).*(rho(index)./atm.r0).^(5./3);
else
out = 2.*(phaseStats.variance(atm)-phaseStats.covariance(rho,atm));
end
end
function out = spectrum(f,atm)
%% SPECTRUM Phase power spectrum density
%
% out = phaseStats.spectrum(f,atm) computes the phase power
% spectrum density from the spatial frequency f and an
% atmosphere object
%
% See also atmosphere
out = (24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6).^2./(2.*pi.^(11./3))).*...
atm.r0.^(-5./3);
out = out.*(f.^2 + 1./atm.L0.^2).^(-11./6);
layers = atm.layer;
out = sum([layers.fractionnalR0]).*out;
end
function out = temporalSpectrum(nu,atm)
%% TEMPORALSPECTRUM Phase temporal power spectrum density
%
% out = phaseStats.temporalSpectrum(nu,atm,zern) computes the
% phase temporal power spectrum density from the temporal
% frequency nu, an atmosphere object and a zernike object
%
% See also atmosphere and zernike
out = zeros(size(nu));
for kLayer = 1:atm.nLayer
atmSlab = slab(atm,kLayer);
[vx,vy] = pol2cart(atmSlab.layer.windDirection,atmSlab.layer.windSpeed);
for k=1:numel(nu)
if vx>eps(atmSlab.layer.windSpeed)
out(k) = out(k) + quadgk( @integrandFy , -Inf, Inf);
else
out(k) = out(k) + quadgk( @integrandFx , -Inf, Inf);
end
end
end
function int = integrandFy(fy)
fx = (nu(k) -fy*vy)/vx;
int = phaseStats.spectrum( hypot(fx,fy) , atmSlab )/vx;
end
function int = integrandFx(fx)
fy = (nu(k) -fx*vx)/vy;
int = phaseStats.spectrum( hypot(fx,fy) , atmSlab )/vy;
end
end
function out = closedLoopVariance(atm,T,tau,gain)
%% SPECTRUM Phase power spectrum density
%
% out = phaseStats.spectrum(f,atm) computes the phase power
% spectrum density from the spatial frequency f and an
% atmosphere object
%
% See also atmosphere
s = @(x) 2*1i*pi*x;
z = @(x) exp(s(x)*T);
G = @(x) ((1-exp(-s(x)*T))./(s(x)*T)).^2.*...
exp(-tau*s(x)).*...
gain./(1-exp(-s(x)*T));
E = @(x) abs(1./(1+G(x)));
figure
nu = logspace(-2,log10(2/T),101);
subplot(1,2,1)
loglog(nu,abs(E(nu)).^2)
xlabel('Hz')
subplot(1,2,2)
loglog(nu,phaseStats.temporalSpectrum(nu,atm),...
nu,phaseStats.temporalSpectrum(nu,atm).*abs(E(nu)).^2)
xlabel('Hz')
drawnow
out = 2*quadgk( @(nu) phaseStats.temporalSpectrum(nu,atm).*abs(E(nu)).^2 , 0 , Inf);
end
function out = symSpectrum(symf)
syms r0 L0
out = (24*gamma(sym(6/5))/5)^(5./6)*...
(gamma(sym(11/6))^2/(2*pi^(11/3)))*r0^(5/3);
out = out*L0^(11/3)*( (symf*L0)^2 + 1 )^(-11/6);
end
function out = otf(rho,atm)
%% OTF Phase optical transfert function
%
% out = phaseStats.otf(rho) compute the phase optical transfert function
% from the baseline rho and an atmosphere object
out = exp(-0.5*phaseStats.structureFunction(rho,atm));
end
function out = psf(f,atm)
%% PSF Phase point spread function
%
% out = phaseStats.psf(f,atm) compute the phase point spread
% from the frequency f and an atmosphere object
fun = @(u) 2.*pi.*quadgk(@(v) psfHankelIntegrandNested(v,u),0,Inf);
out = arrayfun( fun, f);
function y = psfHankelIntegrandNested(x,freq)
y = x.*besselj(0,2.*pi.*x.*freq).*phaseStats.otf(x,atm);
end
end
function out = covarianceMatrix(varargin)
%% COVARIANCEMATRIX Phase covariance matrix
%
% out = phaseStats.covarianceMatrix(rho1,atm) Computes the phase
% auto-covariance matrix from the vector rho1 and an atmosphere
% object
%
% out = phaseStats.covarianceMatrix(rho1,rho2,atm) Computes the phase
% cross-covariance matrix from the vectors rho1 and rho2 and an
% atmosphere object
% Examples:
% covariance matrix on a 1 metre square grid sampled on 16
% pixels :
% [x,y] = meshgrid((0:15)*1/15);
% g = phaseStats.covarianceMatrix(complex(x,y),atm);
% imagesc(g), axis square, colorbar
% covariance matrix on a 1 meter square grid sampled on 16
% pixels with the same grid but displaced of 1 meter:
% [x,y] = meshgrid((0:15)*1/15);
% z = complex(x,y);
% g = phaseStats.covarianceMatrix(z,z+1,atm);
% imagesc(g), axis square, colorbar
%
% See also atmosphere
error(nargchk(2,3,nargin))
rho1 = varargin{1}(:);
if nargin==2
atm = varargin{2};
rho = abs(bsxfun(@minus,rho1,rho1.'));
else
rho2 = varargin{2}(:);
atm = varargin{3};
rho = abs(bsxfun(@minus,rho1,rho2.'));
end
[nRho,mRho] = size(rho);
blockSize = 5000;
if max(nRho,mRho)>blockSize % Memory gentle
l = floor(nRho/blockSize);
le = nRho - l*blockSize;
p = floor(mRho/blockSize);
pe = mRho - p*blockSize;
fprintf(' @(phaseStats.covarianceMatrix)> Memory gentle! (%dX%d blocks)\n',l+1,p+1)
rho = mat2cell( rho, ...
[ones(1,l)*blockSize, le],...
[ones(1,p)*blockSize, pe]);
% out = cell2mat( ...
% cellfun(@(x) phaseStats.covariance(x,atm), rho, 'UniformOutput', false) );
out = cell(size(rho));
for ii=1:numel(rho)
out{ii} = phaseStats.covariance(rho{ii},atm);
end
out = cell2mat(out);
else % Memory intensive
L0r0ratio= (atm.L0./atm.r0).^(5./3);
cst = (24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6)./(2.^(5./6).*pi.^(8./3))).*...
L0r0ratio;
out = ones(size(rho),class(rho)).*(24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6).*gamma(5./6)./(2.*pi.^(8./3))).*L0r0ratio;
index = rho~=0;
u = 2.*pi.*rho(index)./atm.L0;
out(index) = cst.*u.^(5./6).*besselk(5./6,u);
layers = atm.layer;
out = sum([layers.fractionnalR0]).*out;
end
end
function out = covarianceToeplitzMatrix(atm,z1,varargin)
%% COVARIANCETOEPLITZMATRIX Phase covariance matrix
%
% out = phaseStats.covarianceMatrix(atm,rho1) Computes the phase
% auto-covariance matrix from the vector rho1 and an atmosphere
% object
%
% out = phaseStats.covarianceMatrix(atm,rho1,rho2) Computes the phase
% cross-covariance matrix from the vectors rho1 and rho2 and an
% atmosphere object
%
% out = phaseStats.covarianceMatrix(...,'mask',pupil)
%
% Examples:
% covariance matrix on a 1 metre square grid sampled on 16
% pixels :
% [x,y] = meshgrid((0:15)*1/15);
% g = phaseStats.covarianceMatrix(complex(x,y),atm);
% imagesc(g), axis square, colorbar
% covariance matrix on a 1 meter square grid sampled on 16
% pixels with the same grid but displaced of 1 meter:
% [x,y] = meshgrid((0:15)*1/15);
% z = complex(x,y);
% g = phaseStats.covarianceMatrix(z,z+1,atm);
% imagesc(g), axis square, colorbar
%
% See also atmosphere
inputs = inputParser;
inputs.addRequired('atm',@isobject);
inputs.addRequired('z1',@isnumeric);
inputs.addOptional('z2',z1,@isnumeric);
inputs.addParamValue('mask',[],@islogical);
inputs.parse(atm,z1,varargin{:});
atm = inputs.Results.atm;
z1 = inputs.Results.z1;
z2 = inputs.Results.z2;
mask = inputs.Results.mask;
[nz,mz] = size( z1 );
% First Row
r = phaseStats.covariance( abs(z2-z1(1)) , atm );
r = mat2cell( r , nz , ones(mz,1));
% First column in first blocks fow
c = phaseStats.covariance( abs( bsxfun( @minus , z2(1:nz:nz^2), z1(1:nz).' ) ) , atm );
c = mat2cell( c , nz , ones(mz,1) );
% First block rows
rr = cellfun( @(x,y) localToeplitz(x,y) , c , r , 'UniformOutput' , false);
% First Column
c = phaseStats.covariance( abs(z1-z2(1)) , atm );
c = mat2cell( c , nz , ones(mz,1));
% First row in first blocks column
r = phaseStats.covariance( abs( bsxfun( @minus , z1(1:nz:nz^2), z2(1:nz).' ) ) , atm );
r = mat2cell( r , nz , ones(mz,1) );
% First blocks column
cc = cellfun( @(x,y) localToeplitz(x,y) , c , r , 'UniformOutput' , false);
out = cell2mat( localToeplitz(cc,rr) );
out(~mask,:) = [];
out(:,~mask) = [];
function t = localToeplitz(c,r)
% this version works on numeric vector as well as on cells vector
r = r(:); % force column structure
p = length(r);
m = length(c);
x = [r(p:-1:2) ; c(:)]; % build vector of user data
cidx = uint16(0:m-1)';
ridx = uint16(p:-1:1);
subscripts = cidx(:,ones(p,1)) + ridx(ones(m,1),:); % Toeplitz subscripts
t = x(subscripts); % actual data
end
end
function varargout = spatioAngularCovarianceMatrix(sampling,range,atm,srcAC,varargin)
%% SPATIOANGULARCOVARIANCEMATRIX Phase spatio-angular covariance meta matrix
%
% S = spatioAngularCovarianceMatrix(sampling,range,atm,src1)
% computes the spatio-angular auto-correlation meta-matrix of the
% wavefront between all the sources srcAC. The phase is
% sampling with sampling points in the given range and
% propagates through the atmosphere defined by the object atm
%
% C = spatioAngularCovarianceMatrix(sampling,range,atm,src1,src2)
% computes the spatio-angular cross-correlation meta-matrix of the
% wavefront between all src2 and src1. The phase is
% sampling with sampling points in the given range and
% propagates through the atmosphere defined by the object atm
%
% [S,C] = spatioAngularCovarianceMatrix(...) computes both
% auto- and cross-correlation meta-matrix
%
% ... = spatioAngularCovarianceMatrix(...,'mask',mask) restrict
% the sampling within the mask
inputs = inputParser;
inputs.addRequired('sampling',@isnumeric);
inputs.addRequired('range',@isnumeric);
inputs.addRequired('atm',@(x) isa(x,'atmosphere'));
inputs.addRequired('srcAC',@(x) isa(x,'source'));
inputs.addOptional('srcCC',[],@(x) isa(x,'source'));
inputs.addOptional('tipTilt',false,@islogical);
inputs.addParamValue('mask',true(sampling),@islogical);
inputs.addParamValue('lag',0,@isnumeric);
inputs.addParamValue('xyOutput',[],@isnumeric);
inputs.parse(sampling,range,atm,srcAC,varargin{:});
m_srcCC = inputs.Results.srcCC;
tipTilt = inputs.Results.tipTilt;
m_mask = inputs.Results.mask;
m_tau = inputs.Results.lag;
xyOutput= inputs.Results.xyOutput;
[m_x,m_y] = meshgrid( linspace(-1,1,sampling)*range/2 );
m_nGs = length(srcAC);
L0r0ratio= (atm.L0./atm.r0).^(5./3);
m_cst = (24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6)./(2.^(5./6).*pi.^(8./3))).*...
L0r0ratio;
m_cstL0 = (24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6).*gamma(5./6)./(2.*pi.^(8./3))).*L0r0ratio;
m_cstr0 = (24.*gamma(6./5)./5).^(5./6).*...
(gamma(11./6).^2./(2.*pi.^(11./3))).*...
atm.r0.^(-5./3);
m_L0 = atm.L0;
m_nLayer = atm.nLayer;
layers = atm.layer;
m_altitude = [layers.altitude];
m_fr0 = [layers.fractionnalR0];
[m_windVx,m_windVy] = pol2cart([layers.windDirection],[layers.windSpeed]);
m_srcACdirectionVector = cat(2,srcAC.directionVector);
m_srcACheight = [srcAC.height];
if nargout==2
tic
varargout{1} = autoCorrelation(m_x,m_y,m_mask,...
m_nGs,m_srcACdirectionVector,m_srcACheight,...
m_nLayer,m_altitude,m_fr0,...
m_L0,m_cstL0,m_cst);
toc
m_xAC = m_x(m_mask);
m_yAC = m_y(m_mask);
if isempty(xyOutput)
m_xCC = m_xAC;
m_yCC = m_yAC;
else
m_xCC = xyOutput(:,1);
m_yCC = xyOutput(:,2);
end
tic
varargout{2} = crossCorrelation(m_srcCC,m_xAC,m_yAC,m_xCC,m_yCC,...
m_nGs,m_srcACdirectionVector,m_srcACheight,...
m_nLayer,m_altitude,m_fr0,...
m_L0,m_cstL0,m_cst,m_windVx,m_windVy,m_tau);
toc
else
if isempty(m_srcCC)
varargout{1} = autoCorrelation(m_x,m_y,m_mask,...
m_nGs,m_srcACdirectionVector,m_srcACheight,...
m_nLayer,m_altitude,m_fr0,...
m_L0,m_cstL0,m_cst);
else
if tipTilt
varargout{1} = crossCorrelationTT(m_srcCC,m_x,m_y,m_mask,...
m_nGs,m_srcACdirectionVector,m_srcACheight,...
m_nLayer,m_altitude,m_fr0,...
m_L0,m_cstr0,range);
else
m_xAC = m_x(m_mask);
m_yAC = m_y(m_mask);
if isempty(xyOutput)
m_xCC = m_xAC;
m_yCC = m_yAC;
else
m_xCC = xyOutput(:,1);
m_yCC = xyOutput(:,2);
end
varargout{1} = crossCorrelation(m_srcCC,m_xAC,m_yAC,m_xCC,m_yCC,...
m_nGs,m_srcACdirectionVector,m_srcACheight,...
m_nLayer,m_altitude,m_fr0,...
m_L0,m_cstL0,m_cst,m_windVx,m_windVy,m_tau);
end
end
end
function C = crossCorrelationTT(srcTT,x,y,mask,...
nGs,srcACdirectionVector,srcACheight,...
nLayer,altitude,fr0,...
L0,cstr0,D)
fprintf(' -->> Cross-correlation meta-matrix (phase/tip-tilt)!\n')
nSs = length(srcTT);
srcTTdirectionVector = cat(2,srcTT.directionVector);
C = cellfun( @(x) zeros(sum(mask(:))) , cell(nSs,nGs) , 'UniformOutput' , false );
x = x(mask);
y = y(mask);
f02 = 1./L0.^2;
for k=1:nSs*nGs
[kSs,iGs] = ind2sub([nSs,nGs],k);
buf = 0;
for kLayer=1:nLayer
beta = srcACdirectionVector(:,iGs)*altitude(kLayer);
scale = 1 - altitude(kLayer)/srcACheight(iGs);
iZ = complex( x*scale + beta(1) , y*scale + beta(2) );
betaTt = srcTTdirectionVector(:,kSs)*altitude(kLayer);
zTt = iZ.' - complex(betaTt(1),betaTt(2));
rTt = abs(zTt);
oTt = angle(zTt);
Inm = @(r) quadgk( @(f) (f.^2 + f02).^(-11./6).*...
besselj(2,pi.*f.*D).*...
besselj(1,2.*pi.*f.*r),0,Inf);
out = fr0(kLayer)*cstr0.*arrayfun(Inm,rTt).*8/D/scale;
buf = buf + ...
[ out.*cos(oTt) ; out.*sin(oTt) ];
end
C{k} = buf;
end
buf = C;
C = cell(nSs,1);
for k=1:nSs
C{k} = cell2mat(buf(k,:));
end
end
function C = crossCorrelation(srcCC,xAC,yAC,xCC,yCC,...
nGs,srcACdirectionVector,srcACheight,...
nLayer,altitude,fr0,...
L0,cstL0,cst,windVx,windVy,tau)
fprintf(' -->> Cross-correlation meta-matrix!\n')
nSs = length(srcCC);
srcCCdirectionVector = cat(2,srcCC.directionVector);
srcCCheight = [srcCC.height];
C = cellfun( @(x) zeros(length(xCC),length(xAC)) , cell(nSs,nGs) , 'UniformOutput' , false );
cstL0CC = cstL0*ones(length(xCC),length(xAC));
for k=1:nSs*nGs
[kSs,iGs] = ind2sub([nSs,nGs],k);
buf = 0;
for kLayer=1:nLayer
beta = srcACdirectionVector(:,iGs)*altitude(kLayer);
scale = 1 - altitude(kLayer)/srcACheight(iGs);
iZ = complex( xAC*scale + beta(1) , yAC*scale + beta(2) );
betaSs = srcCCdirectionVector(:,kSs)*altitude(kLayer);
scale = 1 - altitude(kLayer)/srcCCheight(kSs);
zSs = complex( ...
xCC*scale + betaSs(1) + windVx(kLayer)*tau, ...
yCC*scale + betaSs(2) + windVy(kLayer)*tau );
rho = abs(bsxfun(@minus,zSs,iZ.'));
out = cstL0CC;
index = rho~=0;
u = 2.*pi.*rho(index)./L0;
out(index) = cst.*u.^(5./6).*besselk(5./6,u);
buf = buf + fr0(kLayer)*out;
end
C{k} = buf;
end
buf = C;
C = cell(nSs,1);
for k=1:nSs
C{k} = cell2mat(buf(k,:));
end
end
function S = autoCorrelation(x,y,mask,...
nGs,srcACdirectionVector,srcACheight,...
nLayer,altitude,fr0,...
L0,cstL0,cst)
fprintf(' -->> Auto-correlation meta-matrix!\n')
kGs = reshape( triu( reshape(1:nGs^2,nGs,nGs) , 1) , 1, []);
kGs(1) = 1;
kGs(kGs==0) = [];
S = cellfun( @(x) zeros(sum(mask(:))) , cell(1,length(kGs)) , 'UniformOutput' , false );
cstL0AC = cstL0*ones(sampling);
for k=1:length(kGs)
[iGs,jGs] = ind2sub( [nGs,nGs] , kGs(k) );
buf = 0;
for kLayer=1:nLayer
% atmSlab = slab(atm,kLayer);
beta = srcACdirectionVector(:,iGs)*altitude(kLayer);
scale = 1 - altitude(kLayer)/srcACheight(iGs);
iZ = complex( x*scale + beta(1) , y*scale + beta(2) );
beta = srcACdirectionVector(:,jGs)*altitude(kLayer);
scale = 1 - altitude(kLayer)/srcACheight(jGs);
jZ = complex( x*scale + beta(1) , y*scale + beta(2) );
z1 = iZ;
z2 = jZ;
[nz,mz] = size( z1 );
% First Row
% r = phaseStats.covariance( abs(z2-z1(1)) , atm );
rho = abs(z2-z1(1));
r = cstL0AC;
index = rho~=0;
u = 2.*pi.*rho(index)./L0;
r(index) = cst.*u.^(5./6).*besselk(5./6,u);
r = mat2cell( fr0(kLayer)*r , nz , ones(mz,1));
% First column in first blocks fow
% c = phaseStats.covariance( abs( bsxfun( @minus , z2(1:nz:nz^2), z1(1:nz).' ) ) , atm );
rho = abs( bsxfun( @minus , z2(1:nz:nz^2), z1(1:nz).' ) );
c = cstL0AC;
index = rho~=0;
u = 2.*pi.*rho(index)./L0;
c(index) = cst.*u.^(5./6).*besselk(5./6,u);
c = mat2cell( fr0(kLayer)*c , nz , ones(mz,1) );
% First block rows
rr = cellfun( @(x,y) myToeplitz(x,y) , c , r , 'UniformOutput' , false);
% First Column
% c = phaseStats.covariance( abs(z1-z2(1)) , atm );
rho = abs(z1-z2(1));
c = cstL0AC;
index = rho~=0;
u = 2.*pi.*rho(index)./L0;
c(index) = cst.*u.^(5./6).*besselk(5./6,u);
c = mat2cell( fr0(kLayer)*c , nz , ones(mz,1));
% First row in first blocks column
% r = phaseStats.covariance( abs( bsxfun( @minus , z1(1:nz:nz^2), z2(1:nz).' ) ) , atm );
rho = abs( bsxfun( @minus , z1(1:nz:nz^2), z2(1:nz).' ) );
r = cstL0AC;
index = rho~=0;
u = 2.*pi.*rho(index)./L0;
r(index) = cst.*u.^(5./6).*besselk(5./6,u);
r = mat2cell( fr0(kLayer)*r , nz , ones(mz,1) );
% First blocks column
cc = cellfun( @(x,y) myToeplitz(x,y) , c , r , 'UniformOutput' , false);
out = cell2mat( myToeplitz(cc,rr) );
out(~mask,:) = [];
out(:,~mask) = [];
buf = buf + out;
end
S{k} = buf;
end
buf = S;
S = cellfun( @(x) zeros(sum(mask(:))) , cell(nGs) , 'UniformOutput' , false );
S(kGs) = buf;
S(1:nGs+1:nGs^2) = S(1,1);
S = cell2mat(S);
S = triu(S,1)+triu(S)';
end
end
function out = upwardJitter(L,D,atm,zSrc)
%% UPWARDJITTER Upward propagated beam motion
%
% varJit = upwardJitter(L,D,atm) computes the variance of the
% upward propagated beam motion at a target at distance L from
% source, the beam is lauched from a telescope of diameter D
% and is propagated through the atmosphere object atm
%
% varJit = upwardJitter(L,D,atm,zSrc) computes the variance of
% the upward propagated beam motion at a target at distance L
% from source at zSrc, the beam is lauched from a telescope of
% diameter D and is propagated through the atmosphere object
% atm
if nargin<4
zSrc = L;
end
k0 = 2*pi/atm.wavelength;
fun = @(f,z,atmSlab) f.*phaseStats.spectrum(f,atmSlab).*...
(L-z).^2.*...
(16/(k0.*(1-z./zSrc).*D)).^2.*...
(besselj(2,2*pi*f*(1-z./zSrc)*D/2)./(2*pi*f*(1-z./zSrc)*D/2)).^2;
out = 0;
for kLayer=1:atm.nLayer
atmSlab = slab(atm,kLayer);
z = atmSlab.layer.altitude;
out = out + 2.*pi.*quadgk( @(x) fun(x,z,atmSlab) , 0 , Inf);
end
end
function out = zernikeVariance(zern,atm)
%% ZERNIKEVARIANCE Zernike coefficients variance
%
% out = variance(modes,atmosphere) computes the
% variance of Zernike coefficients from the modes and the
% atmosphere object
%
% out = variance(zernike,atmosphere) computes the
% variance of Zernike coefficients from the Zernike polynomials
% object and the atmosphere object
%
% Example:
% atm = atmosphere(photometry.V,0.15,30);
% tel = telescope(10);
% modes = 1:15;
% figure
% semilogy(modes,phaseStats.zernikeVariance(modes,atm),'--.')
% xlabel('Zernike modes')
% ylabel('Variance [rd^2]')
%
% See also zernike, atmosphere
if ~isa(zern,'zernike')
zern = zernike(zern);
end
r0 = atm.r0;
L0 = atm.L0;
D = zern.D;
jv = zern.j;
nv = zern.n;
mv = zern.m;
nv0 = nv;
index = diff(nv)~=0;
jv = [jv(index) jv(end)];
mv = [mv(index) mv(end)];
nv = [nv(index) nv(end)];
nf = length(nv);
out = zeros(length(zern.j),1);
for cpt = 1:nf
j = jv(cpt);
n = nv(cpt);
m = mv(cpt);
out(nv0==n,1) = zernCovCoef(r0,L0,D,j,j,n,m,n,m);
end
function out = zernCovCoef(r0,L0,D,i,j,ni,mi,nj,mj)
if (mi==mj) && (rem(abs(i-j),2)==0 || ((mi==0) && (mj==0)))
if L0==Inf
if i==1 && j==1
out = Inf;
else
out = (gamma(11./6).^2.*gamma(14./3)./(2.^(8./3).*pi)).*(24.*gamma(6./5)./5).^(5./6).*...
(D./r0).^(5./3).*sqrt((ni+1).*(nj+1)).*(-1).^((ni+nj-mi-mj)./2).*...
newGamma(-5./6+(ni+nj)./2,...
[23./6+(ni+nj)./2 17./6+(ni-nj)./2 17./6+(nj-ni)./2]);
end
else
out = (4.*gamma(11./6).^2./pi.^(14./3)).*(24.*gamma(6./5)./5).^(5./6).*...
(L0./r0).^(5./3).*(L0./D).^2.*...
sqrt((ni+1).*(nj+1)).*(-1).^((ni+nj-mi-mj)./2).*...
UnParamEx4q2(0,ni+1,nj+1,11./6,pi.*D./L0);
end
else
out = 0;
end
function out = newGamma(a,b)
% NEWGAMMA Computes the function defined by Eq.(1.18) in R.J. Sasiela's book :
% Electromagnetic Wave Propagation in Turbulence, Springer-Verlag.
% out = newGamma(a,b)
out = prod(gamma(a))./prod(gamma(b));
end
function out = UnParamEx4q2(mu,alpha,beta,p,a)
% UNPARAMEX4Q2 Computes the integral given by the Eq.(2.33) of the thesis
% of R. Conan (Modelisation des effets de l'echelle externe de coherence
% spatiale du front d'onde pour l'Observation a Haute Resolution Angulaire
% en Astronomie, University of Nice-Sophia Antipolis, October 2000)
% http://www-astro.unice.fr/GSM/Bibliography.html#thesis
a1 = [(alpha+beta+1)./2 (2+mu+alpha+beta)./2 (mu+alpha+beta)./2];
b1 = [1+alpha+beta 1+alpha 1+beta];
a2 = [(1-mu)./2+p 1+p p];
b2 = [1+(alpha+beta-mu)./2+p 1+(alpha-beta-mu)./2+p 1+(beta-alpha-mu)./2+p];
out = (1./(2.*sqrt(pi).*gamma(p))).*(...
newGamma([a1 p-(mu+alpha+beta)./2],b1).*a.^(mu+alpha+beta).*...
pochammerSeries(3,5,a1,[1-p+(mu+alpha+beta)./2 b1 1],a.^2) + ...
newGamma([(mu+alpha+beta)./2-p a2],b2).*a.^(2.*p).*...
pochammerSeries(3,5,a2,[1-(mu+alpha+beta)./2+p b2 1],a.^2));
function out = pochammerSeries(p,q,a,b,z,tol,nmax)
% POCHAMMERSERIES Computes power series in Pochammer notation
% pochammerSeries(p,q,a,b,z)
% pochammerSeries(p,q,a,b,z,tol)
% pochammerSeries(p,q,a,b,z,[],nmax)
% pochammerSeries(p,q,a,b,z,tol,nmax)
if (p==(q+1) && abs(z)<1) || (abs(z)==1 && real(sum(a)-sum(b))<0) || p<(q+1)
if p==length(a) && q==length(b)
switch nargin
case 6
nmax = 1e3;
case 7
if isempty(tol)
tol = 1e-6;
end
otherwise
tol = 1e-6;
nmax = 1e3;
end
out = zeros(size(z));
indz = find(z==0);
if ~isempty(indz)
out(indz) = 1;
end
indnz = find(z~=0);
if ~isempty(indnz)
z = z(indnz);
ck = 1;
step = Inf;
k = 0;
som = ck;
while (k<=nmax) && (step>tol)
ckp1 = prod(a+k).*z.*ck./prod(b+k);
step = abs(abs(ck)-abs(ckp1));
som = som + ckp1;
k = k+1;
ck = ckp1;
end
if step>tol
warning('pochammerSeries','Maximum iteration reached before convergence')
end
out(indnz) = som;
end
else
error('p and q must be the same length than vectors a and b, respectively')
end
else
error('This generalized hypergeometric function doesn''t converge')
end
end
end
end
end
function out = zernikeCovariance(zern,atm)
%% ZERNIKECOVARIANCE Zernike coefficients covariance
%
% out = phaseStats.zernikeCovariance(modes,atmosphere)
% computes the covariance matrix of Zernike coefficients from
% the modes and the atmosphere object
%
% out = phaseStats.zernikeCovariance(zernike,atmosphere) computes the
% covariance matrix of Zernike coefficients from the Zernike
% polynomials object and the atmosphere object
%
% Example:
% atm = atmosphere(photometry.V,0.15,30);
% tel = telescope(10);
% modes = 1:15;
% figure
% spy(phaseStats.zernikeCovariance(modes,atm,tel))
%
% See also zernike, atmosphere
if ~isa(zern,'zernike')
zern = zernike(zern);
end
r0 = atm.r0;
L0 = atm.L0;
D = zern.D;
% -- covariance --
n = zern.n;
m = zern.m;
[i,j] = meshgrid(zern.j);
[ni,nj] = meshgrid(n);
[mi,mj] = meshgrid(m);
index = (mi==mj) & (rem(abs(i-j),2)==0 | ((mi==0) & (mj==0)));
index = logical(index - diag(diag(index)));
i = i(index);
j = j(index);