-
Notifications
You must be signed in to change notification settings - Fork 11
/
Copy pathPseudoCl.F90
1833 lines (1512 loc) · 75.2 KB
/
PseudoCl.F90
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
!Functions for calculating coupling matrices and covariances, etc.
!Desgined for use on a cluster 2+ CPUs using MPI, not so well tested on single CPU
!This version AL: April 2008 - Dec 2010
module PseudoCl
use healpix_types, ONLY: SP,DP,I4B,SPC
use HealpixObj
use MpiStuff
USE spinalm_tools
implicit none
Type TBeam
logical :: beam_transfer
real(dp) :: fwhm
real(dp), pointer :: beam(:) => NULL()
end type TBeam
Type TCrossBeamSet
Type(TBeam), allocatable :: Beams(:,:)
end Type
Type TChannel
character(LEN=16) :: Name
character(LEN=256) :: noise_file
integer :: Count !Number of detectors
Type (TStringList) :: DetectorNames
logical, pointer :: detector_has_pol(:) => NULL() !true if TQU, false if just T
real(dp), pointer :: sig0(:) => NULL()
real(dp) :: Ghz
real(dp) :: PtSrcA !non-beam-smoothed C_l of fiducial point sources
real(dp) :: ENoiseFac
Type(TBeam), pointer :: DetectorBeams(:) => NULL()
Type(HealpixMap), pointer :: DetectorYearNoiseMaps(:,:) => NULL()
Type(TBeam) :: Beam
Type(HealpixMap) :: NoiseMap
Type(HealpixPower), pointer :: NoiseP(:) => NULL()
end Type TChannel
Type TCouplingMatrix
integer lmin
integer lmax
integer nl
logical has_pol
real(dp), dimension(:,:), pointer :: T => NULL() , X => NULL() , &
EE => NULL() , EB => NULL()
real(dp), dimension(:,:), pointer :: InvT => NULL() ,InvX => NULL() , &
InvEE => NULL() , InvEB => NULL()
end type TCouplingMatrix
Type TCovMat
real(dp), dimension(:,:), pointer :: C => NULL()
end Type TCovMat
Type TCovMatArray
integer ncl
Type(TCovMat), dimension(:,:), pointer :: Cov => NULL()
end Type TCovMatArray
Type TCovMatPolArray
integer vec_size
Type(TCovMatArray), dimension(:,:), pointer :: Pol => NULL()
end Type TCovMatPolArray
Type TCovMatSet
integer n, lmin, lmax
Type(TCovMat), dimension(:), pointer :: Cov => NULL()
end Type TCovMatSet
Type TClArray
Type(HealpixPower), dimension(:), pointer :: P => NULL()
end Type TClArray
Type TCouplingMatrixArray
Type(TCouplingMatrix), dimension(:), pointer :: M => NULL()
end Type TCouplingMatrixArray
Real(dp):: NoiseScale
!Avoid very small numbers is noise*w_i*w_j maps;
!set by WeightsToCovPowers
!used by GetFullCovariance
real(dp) :: cross_noise_scale = 1.d0
contains
subroutine TBeam_Free(B)
Type(TBeam) :: B
integer status
deallocate(B%Beam,stat = status)
nullify(B%Beam)
end subroutine TBeam_Free
subroutine TBeam_Assign(B, Bin)
Type(TBeam) B, Bin
call TBeam_Free(B)
B%beam_transfer = Bin%beam_transfer
B%fwhm = Bin%fwhm
if (B%beam_transfer) then
allocate(B%beam(0:size(Bin%beam)-1))
B%beam = Bin%Beam
end if
end subroutine TBeam_Assign
subroutine TBeam_SetGaussian(B, lmax)
Type(TBeam) :: B
integer, intent(in) :: lmax
real xlc, sigma2
integer l
B%beam_transfer = .true.
allocate(B%beam(0:lmax))
B%beam = 1
xlc= 180*sqrt(8.*log(2.))/HO_pi
sigma2 = (B%fwhm/xlc)**2
do l=1,lmax
B%beam(l) = exp(-l*(l+1)*sigma2/2)
end do
end subroutine TBeam_SetGaussian
subroutine TBeam_ReadFile(B, filename, lmax, col)
Type(TBeam) :: B
character(LEN=*), intent(in) :: filename
integer, intent(in), optional :: col
integer, intent(in) :: lmax
integer i,l, beamcol
real(dp) inB
character(LEN=1024*8) :: inline
real(dp), allocatable :: cols(:)
if (present(col)) then
beamcol=col
else
beamcol=2
end if
allocate(cols(col-1))
B%beam_transfer = .true.
allocate(B%beam(0:lmax))
call OpenTxtFile(filename,1)
do
read (1,'(a)', err = 500, end =500) inline
if (inline(1:1) /= '#' .and. inline/='') then
read(inline,*) l, cols
if (l<=lmax) B%beam(l) = cols(beamcol-1)
end if
end do
500 close(1)
if (l<lmax) call MpiStop('TBeam_ReadFile: not high enough l in file')
end subroutine TBeam_ReadFile
subroutine TBeam_PowerSmooth(B, P, sgn)
Type(TBeam) :: B
Type (HealpixPower)::P
integer, intent(in) :: sgn
if (P%lmax > size(B%Beam)-1) call MpiStop('TBeam_PowerSmooth: lmax mismatch')
if (B%beam_transfer) then
call HealpixPower_Smooth_Beam(P,B%beam,sgn)
else
call HealpixPower_Smooth(P,B%fwhm,sgn)
end if
end subroutine TBeam_PowerSmooth
subroutine TBeam_PowerSmooth2(B1,B2, P, sgn)
Type(TBeam) :: B1, B2
Type (HealpixPower)::P
integer, intent(in) :: sgn
if (P%lmax > size(B1%Beam)-1) call MpiStop('TBeam_PowerSmooth2: lmax mismatch')
if (.not. B1%beam_transfer) call MpiStop('TBeam_PowerSmooth2: want transfers now')
call HealpixPower_Smooth_Beam2(P,B1%beam,B2%beam,sgn)
end subroutine TBeam_PowerSmooth2
subroutine TCouplingMatrix_Free(M)
Type(TCouplingMatrix) :: M
integer status
deallocate(M%T,stat=status)
deallocate(M%X, M%EE, M%EB, stat=status)
deallocate(M%InvT,stat=status)
deallocate(M%InvX, M%InvEE,M%InvEB, stat=status)
call TCouplingMatrix_Nullify(M)
end subroutine TCouplingMatrix_Free
subroutine TCouplingMatrix_ArrayFree(C)
Type(TCouplingMatrix) :: C(:)
integer i
do i=1,size(C)
call TCouplingMatrix_Free(C(i))
end do
end subroutine TCouplingMatrix_ArrayFree
subroutine TCouplingMatrix_Nullify(M)
Type(TCouplingMatrix) :: M
nullify(M%T,M%X,M%EE,M%EB)
nullify(M%InvT,M%InvX,M%InvEE,M%InvEB)
end subroutine TCouplingMatrix_Nullify
subroutine PseudoCl_GetCouplingInversesArr(M, n)
!Get M^{-1} by inverting symmetric form and putting in l factors
use MatrixUtils
integer, intent(in) :: n
Type(TCouplingMatrix) :: M(n)
Type(TCouplingMatrix), pointer :: AM
integer lmin,lmax,l
real(dm), dimension(:,:), allocatable :: tmp,Inv11
integer MpiID, MpiSize, i,j
integer params(3)
logical haspol
Type(TMatrixType), dimension(:), allocatable :: Arr
Type(TMatrixType) dummyArr(1)
integer nmat
call GetMpiStat(MpiID, MpiSize)
if (MpiID==0) then
nmat = n
lmin=M(1)%lmin
lmax=M(1)%lmax
if (M(1)%has_pol) nmat=n*3
print *,'nmat = ',nmat
allocate(Arr(nmat))
do i=1,n
allocate(Arr(i)%M(lmin:lmax,lmin:lmax))
Arr(i)%M = M(i)%T
if (M(i)%has_pol) then
print *,'inverting pol'
allocate(Arr(i+n)%M(lmin:lmax,lmin:lmax))
Arr(i+n)%M = M(i)%EE
allocate(Arr(i+2*n)%M(lmin:lmax,lmin:lmax))
Arr(i+2*n)%M = M(i)%X
end if
end do
call Matrix_InverseArrayMPI(Arr,nmat)
do i=1,n
M(i)%InvT => Arr(i)%M
do l=lmin,lmax
M(i)%InvT(l,:)=M(i)%InvT(l,:)/(2*l+1)
end do
if (M(i)%has_pol) then
M(i)%InvEE => Arr(i+n)%M
M(i)%InvX => Arr(i+2*n)%M
allocate(inv11(lmin:lmax,lmin:lmax))
allocate(tmp(lmin:lmax,lmin:lmax))
Inv11= M(i)%InvEE
!M%InvEE = M%EE - matmul(M%EB,matmul(Inv11,M%EB))
call Matrix_Mult_SymmRight(Inv11,M(i)%EB,tmp)
call Matrix_Mult(M(i)%EB,tmp,M(i)%InvEE)
do l=lmin,lmax
M(i)%InvEE(:,l) = M(i)%EE(:,l) - M(i)%InvEE(:,l)
end do
call Matrix_inverse(M(i)%InvEE)
print *,'getting InvEB'
call Matrix_Mult_SymmRight(M(i)%InvEE,M(i)%EB,tmp)
allocate(M(i)%InvEB(lmin:lmax,lmin:lmax))
call Matrix_Mult_SymmRight(tmp,inv11,M(i)%InvEB,-1._dm)
deallocate(Inv11,tmp,Arr)
do l=lmin,lmax
M(i)%InvX(l,:)=M(i)%InvX(l,:)/(2*l+1)
M(i)%InvEE(l,:)=M(i)%InvEE(l,:)/(2*l+1)
M(i)%InvEB(l,:)=M(i)%InvEB(l,:)/(2*l+1)
end do
end if
end do
else
call Matrix_InverseArrayMPI(DummyArr,nmat)
end if
end subroutine PseudoCl_GetCouplingInversesArr
subroutine TCovMatArray_Free(C)
Type(TCovmatArray) :: C
integer i,j, err
do i=1, C%ncl
do j=1, C%ncl
deallocate(C%Cov(i,j)%C, stat=err)
end do
end do
deallocate(C%Cov)
end subroutine TCovMatArray_Free
subroutine TCouplingMatrix_Write(M, fname,i)
use AMLUtils
Type(TCouplingMatrix) :: M
character(LEN=*), intent(in) :: fname
integer, intent(in) :: i
logical B1,B2
call CreateFile(fname, i ,'UNFORMATTED')
B1 = associated(M%T)
B2 = associated(M%InvT)
write (i) M%lmin,M%lmax,M%nl,M%has_pol,B1,B2
if (B1) then
write(i) M%T
if (M%has_pol) then
write(i) M%X
write(i) M%EE
write(i) M%EB
end if
end if
if (B2) then
write(i) M%InvT
if (M%has_pol) then
write(i) M%InvX
write(i) M%InvEE
write(i) M%InvEB
end if
end if
Close(i)
end subroutine TCouplingMatrix_Write
subroutine TCouplingMatrix_Read(M, fname,i)
use AMLUtils
Type(TCouplingMatrix) :: M
character(LEN=*), intent(in) :: fname
integer, intent(in) :: i
integer lmin,lmax
logical B1,B2
call TCouplingMatrix_Free(M)
call OpenFile(fname, i ,'UNFORMATTED')
read (i) lmin,lmax,M%nl,M%has_pol,B1,B2
M%lmin=lmin
M%lmax=lmax
if (B1) then
allocate(M%T(lmin:lmax,lmin:lmax))
read(i) M%T
if (M%has_pol) then
allocate(M%X(lmin:lmax,lmin:lmax))
allocate(M%EE(lmin:lmax,lmin:lmax))
allocate(M%EB(lmin:lmax,lmin:lmax))
read(i) M%X
read(i) M%EE
read(i) M%EB
end if
else
nullify(M%T,M%X,M%EE,M%EB)
end if
if (B2) then
allocate(M%InvT(lmin:lmax,lmin:lmax))
read(i) M%InvT
if (M%has_pol) then
allocate(M%InvX(lmin:lmax,lmin:lmax))
allocate(M%InvEE(lmin:lmax,lmin:lmax))
allocate(M%InvEB(lmin:lmax,lmin:lmax))
read(i) M%InvX
read(i) M%InvEE
read(i) M%InvEB
end if
else
nullify(M%InvT,M%InvEE,M%InvEB,M%InvT)
end if
Close(i)
end subroutine TCouplingMatrix_Read
subroutine SetArrayInterlacedSym(MMixed,Mout, nthreads)
real(dp) :: Mout(:,:)
real(dp), intent(in) :: MMixed(:,:,:)
integer, intent(in) :: nthreads
integer id, i, j,colix
do id = 0, nthreads-1
colix=0
do i = 1+id, size(Mout,DIM=2), nthreads
colix=colix+1
Mout(1:i,i) = MMixed(1:i,colix,id+1)
end do
end do
do i = 1, size(Mout,DIM=2)
do j=1, i-1
Mout(i,j) = Mout(j,i)
end do
end do
end subroutine SetArrayInterlacedSym
subroutine AssignSym(Mout, Minp)
real(dp) :: Mout(:,:), Minp(:,:)
integer i,j
do i = 1, size(Mout,DIM=2)
do j=1, i
Mout(i,j) = Minp(j,i)
Mout(j,i) = Minp(j,i)
end do
end do
end subroutine AssignSym
subroutine PseudoCl_GetCouplingMatrixArr(M_in, P, inlmin, inlmax, indopol, inncl, inpolweights, inpolspin)
!Get coupling matrix from power spectrum of the weight function
!For the moment assume weight same for all maps
!Returns the symmetric form M_{l1 l2}/(2*l2+1)
use AMLUtils
use spinalm_tools
Type(TCouplingMatrix),dimension(:), target :: M_in
Type(TCouplingMatrix),dimension(:), pointer :: M
Type(HealpixPower),dimension(:,:) :: P
logical, intent(in) :: indopol
integer, intent(in) :: inlmin,inlmax, inncl
logical, intent(in) :: inpolweights
integer, intent(in), optional :: inpolspin
integer ncl,lmax, npolweights
logical dopol
real(dp), dimension(:,:,:), allocatable :: W
real(dp), dimension(:,:), allocatable :: tmp
real(dp), dimension(:,:,:), allocatable :: blocks
real(dp), dimension(:), allocatable :: threejj0, threejj2, sthreejj0, sthreejj2,sthreejj20
integer :: lmin, nl, l1, l2, lplus, lminus, Plmax
integer clix
integer MpiID, MpiSize
integer colix, ierr, ncols, polspin
integer id, params(6), TPix, PPix
double precision IniTime
!Assume loads of memory
IniTime = GeteTime()
call GetMpiStat(MpiID, MpiSize)
params(1) = inlmax
params(2) = inncl
if (indopol) then
params(3)=1
else
params(3)=0
end if
params(4) = inlmin
if (inpolweights) then
params(5) = 3
else
params(5)= 1
end if
if (.not. present(inpolspin)) then
params(6)=2
else
params(6)=inpolspin
end if
#ifdef MPIPIX
call MPI_BCAST(params,SIze(params),MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
#endif
lmax=params(1)
ncl=params(2)
dopol = params(3)==1
lmin=params(4)
npolweights = params(5)
polspin = params(6)
if (npolweights==1) then
TPix=1
PPix=1
else
TPix=2
PPix=3
end if
if (MpiID>0) then
allocate(M(ncl))
else
M => M_in
end if
allocate(W(0:2*lmax, ncl, npolweights))
W=0
allocate(threejj0(0:2*lmax), sthreejj0(0:2*lmax))
if (dopol) then
allocate(threejj2(0:2*lmax), sthreejj2(0:2*lmax), sthreejj20(0:2*lmax))
end if
ncols = (lmax - lmin + MpiSize) / MpiSize !Max number
do clix = 1, ncl
if (MpiID==0) then
! print *,'free'
call TCouplingMatrix_Nullify(M(clix)) !Currently seg faults if try to free, usual intel bug
if (lmax > P(clix,1)%lmax) call MpiStop('PseudoCl_GetCouplingMatrix: lmax in power spectrum too low')
if (lmax*2 > P(clix,1)%lmax) write(*,*) 'Warning: PseudoCl_GetCouplingMatrix really needs 2*lmax '
end if
M(clix)%has_pol = dopol
M(clix)%lmin = lmin
M(clix)%lmax = lmax
M(clix)%nl = lmax-lmin +1
allocate(M(clix)%T(lmin:lmax,ncols))
M(clix)%T=0
if (dopol) then
allocate(M(clix)%X(lmin:lmax,ncols))
M(clix)%X=0
allocate(M(clix)%EE(lmin:lmax,ncols))
M(clix)%EE=0
allocate(M(clix)%EB(lmin:lmax,ncols))
M(clix)%EB=0
end if
if (MpiID==0) then
do l1 = 0, min(P(clix,1)%lmax,2*lmax)
W(l1,clix,1) = (2*l1+1)*P(clix,1)%Cl(l1,C_T)/(4*pi)
if (npolweights>1) then
W(l1,clix,2) = (2*l1+1)*P(clix,2)%Cl(l1,C_T)/(4*pi)
W(l1,clix,3) = (2*l1+1)*P(clix,3)%Cl(l1,C_T)/(4*pi)
end if
end do
end if
end do
#ifdef MPIPIX
call MPI_BCAST(W,SIze(W),MPI_DOUBLE_PRECISION, 0, MPI_COMM_WORLD, ierr)
#endif
colix =0
do l1 = lmin + MpiID, lmax, MpiSize
colix = colix+1
if (colix > ncols) stop 'colix > ncols'
! print *,MpiID, 'colix',colix
do l2 = lmin, l1
lplus = l1+l2
lminus = abs(l1-l2)
call GetThreeJs(threejj0(lminus:),l1,l2,0,0)
sthreejj0(lminus:lplus) = threejj0(lminus:lplus)**2
if (dopol) then
!note that lminus is correct, want max(abs(l1-l2),abs(m1)) where m1=0 here
call GetThreeJs(threejj2(lminus:),l1,l2,-polspin,polspin)
sthreejj2(lminus:lplus) = threejj2(lminus:lplus)**2
sthreejj20(lminus:lplus:2) = threejj0(lminus:lplus:2)*threejj2(lminus:lplus:2)
end if
do clix = 1, ncl
M(clix)%T(l2,colix) = sum(W(lminus:lplus,clix,1)*sthreejj0(lminus:lplus))
if (dopol) then
M(clix)%X(l2,colix) = sum(W(lminus:lplus:2,clix,TPix)*sthreejj20(lminus:lplus:2))
M(clix)%EE(l2,colix) = sum(W(lminus:lplus:2,clix,PPix)*sthreejj2(lminus:lplus:2))
M(clix)%EB(l2,colix) = sum(W(lminus+1:lplus:2,clix,PPix)*sthreejj2(lminus+1:lplus:2))
end if
end do
end do
end do
deallocate(W)
deallocate(threejj0,sthreejj0)
if (dopol) deallocate(threejj2,sthreejj2,sthreejj20)
#ifdef MPIPIX
allocate(Blocks(lmin:lmax, ncols,MpiSize))
#else
allocate(tmp(lmin:lmax,lmin:lmax))
#endif
do clix = 1, ncl
#ifdef MPIPIX
call MPI_GATHER(M(clix)%T,size(M(clix)%T),MPI_DOUBLE_PRECISION,blocks,size(M(clix)%T),MPI_DOUBLE_PRECISION,0,MPI_COMM_WORLD,ierr)
deallocate(M(clix)%T)
if (MpiID==0) then
allocate(M(clix)%T(lmin:lmax,lmin:lmax))
call SetArrayInterlacedSym(blocks, M(clix)%T,MpiSize)
end if
if (dopol) then
call MPI_GATHER(M(clix)%X,size(M(clix)%X),MPI_DOUBLE_PRECISION,Blocks,size(M(clix)%X),MPI_DOUBLE_PRECISION,0,MPI_COMM_WORLD,ierr)
deallocate(M(clix)%X)
if (MpiID==0) then
allocate(M(clix)%X(lmin:lmax,lmin:lmax))
call SetArrayInterlacedSym(blocks, M(clix)%X,MpiSize)
end if
call MPI_GATHER(M(clix)%EE,size(M(clix)%EE),MPI_DOUBLE_PRECISION,Blocks,size(M(clix)%EE),MPI_DOUBLE_PRECISION,0,MPI_COMM_WORLD,ierr)
deallocate(M(clix)%EE)
if (MpiID==0) then
allocate(M(clix)%EE(lmin:lmax,lmin:lmax))
call SetArrayInterlacedSym(Blocks, M(clix)%EE,MpiSize)
end if
call MPI_GATHER(M(clix)%EB,size(M(clix)%EB),MPI_DOUBLE_PRECISION,Blocks,size(M(clix)%EB),MPI_DOUBLE_PRECISION,0,MPI_COMM_WORLD,ierr)
deallocate(M(clix)%EB)
if (MpiID==0) then
allocate(M(clix)%EB(lmin:lmax,lmin:lmax))
call SetArrayInterlacedSym(Blocks, M(clix)%EB,MpiSize)
end if
end if
#else
tmp = M(clix)%T
deallocate(M(clix)%T)
allocate(M(clix)%T(lmin:lmax,lmin:Lmax))
call AssignSym(M(clix)%T,tmp)
! M(clix)%T = tmp
if (dopol) then
tmp = M(clix)%X
deallocate(M(clix)%X)
allocate(M(clix)%X(lmin:lmax,lmin:Lmax))
call AssignSym(M(clix)%X,tmp)
tmp = M(clix)%EE
deallocate(M(clix)%EE)
allocate(M(clix)%EE(lmin:lmax,lmin:Lmax))
call AssignSym(M(clix)%EE, tmp)
tmp = M(clix)%EB
deallocate(M(clix)%EB)
allocate(M(clix)%EB(lmin:lmax,lmin:Lmax))
call AssignSym(M(clix)%EB, tmp)
end if
#endif
end do
#ifdef MPIPIX
if (MpiID>0) deallocate(M)
deallocate(Blocks)
#else
deallocate(tmp)
#endif
if (MpiID==0) print *,'coupling time:', GeteTime() - IniTime
end subroutine PseudoCl_GetCouplingMatrixArr
subroutine PseudoCl_GetCHat(M, PCls, Phat)
!We only get TE where map index of E is >= that of T
use MatrixUtils
integer lmax
Type(TCouplingMatrix) :: M
Type(HealpixAllCl) :: Pcls
Type(HealpixPower) :: Phat
real(dm) tmp1(M%lmin:M%lmax),tmp2(M%lmin:M%lmax),tmpB(M%lmin:M%lmax)
if (.not. associated(M%InvT)) then
Stop 'PseudoCl_GetCHat: need inverses'
end if
call HealpixPower_Init(Phat,M%lmax, M%has_pol .and. size(PCls%Cl,2)+size(PCls%Cl,3)>2)
tmp1=PCls%Cl(M%lmin:M%lmax,1,1)
call Matrix_MulVec(M%InvT,tmp1,tmp2)
Phat%Cl(M%lmin:M%lmax,C_T) = tmp2
if (M%has_pol) then
if (size(PCls%Cl,2)>1) then
tmp1=PCls%Cl(M%lmin:M%lmax,2,1)
call Matrix_MulVec(M%InvX,tmp1,tmp2)
Phat%Cl(M%lmin:M%lmax,C_C)=tmp2
if (size(PCls%Cl,3)>1) then
tmp1 = PCls%Cl(M%lmin:M%lmax,2,2)
tmpB = PCls%Cl(M%lmin:M%lmax,3,3)
call Matrix_MulVec(M%InvEE,tmp1,tmp2)
call Matrix_MulVec(M%InvEB,tmpB,tmp2,1._dm,1._dm)
Phat%Cl(M%lmin:M%lmax,C_E)=tmp2
call Matrix_MulVec(M%InvEE,tmpB,tmp2)
call Matrix_MulVec(M%InvEB,tmp1,tmp2,1._dm,1._dm)
Phat%Cl(M%lmin:M%lmax,C_B)=tmp2
end if
end if
end if
end subroutine PseudoCl_GetCHat
subroutine PseudoCl_GetSimpleCovariance(M, Cov, PFid, vec_size, NoiseP, fwhm)
!Simple 1/0 mask with noise added into the theory Cl
!Not used
use MatrixUtils
Type(TCouplingMatrix) :: M
Type(HealpixPower) :: PFid, NoiseP
Type(TCovMat) :: Cov
integer, intent(in) :: vec_size
real(dp), intent(in) :: fwhm
integer nl,i, status
real fac
real rootT(M%lmin:M%lmax), rootE(M%lmin:M%lmax), rootB(M%lmin:M%lmax)
real(dp), dimension(:,:), allocatable :: Minv,tmp
nl = M%lmax-M%lmin+1
if (.not. associated(M%InvT)) then
call MpiStop('Need inverses')
end if
PFid%Cl = PFid%Cl + NoiseP%Cl
deallocate(Cov%C,stat=status)
allocate(Cov%C(nl*vec_size,nl*vec_size))
Cov%C = 0
allocate(Minv(nl*vec_size,nl*vec_size))
rootT = sqrt(PFid%Cl(M%lmin:M%lmax,C_T))
rootE = sqrt(PFid%Cl(M%lmin:M%lmax,C_E))
rootB = sqrt(PFid%Cl(M%lmin:M%lmax,C_B))
print *,'doing T'
do i=1,nl
Cov%C(1:nl,i) = 2*PFid%Cl(M%lmin+i-1,C_T)*M%T(:,M%lmin+i-1)*PFid%Cl(M%lmin:M%lmax,C_T)
end do
print *,'doing Minv'
Minv(1:nl,1:nl) = M%invT
print *,'doing pol'
if (M%has_pol) then
!ordering is T X E B
if (vec_size /= 3 .and. vec_size /= 4) &
stop 'PseudoCl_GetCovariance: vec_size must be 3 (TT, TE, EE) or 4 (+BB)'
do i=1,nl
!XX
fac = rootE(M%lmin+i-1)*rootT(M%lmin+i-1)
Cov%C(nl+1:2*nl,nl+i) = fac*M%X(:,M%lmin+i-1)*rootE(M%lmin:M%lmax)*rootT(M%lmin:M%lmax) &
+ 2*PFid%Cl(M%lmin+i-1,C_C)*M%T(:,M%lmin+i-1)*PFid%Cl(M%lmin:M%lmax,C_C)
!EE
Cov%C(2*nl+1:3*nl,2*nl+i) = 2*PFid%Cl(M%lmin+i-1,C_E)*M%EE(:,M%lmin+i-1)*PFid%Cl(M%lmin:M%lmax,C_E) &
+2*PFid%Cl(M%lmin+i-1,C_B)*M%EB(:,M%lmin+i-1)*PFid%Cl(M%lmin:M%lmax,C_B)
!TX
Cov%C(1:nl,nl+i) = rootT(M%lmin+i-1)*rootT(M%lmin:M%lmax)*&
(PFid%Cl(M%lmin:M%lmax,C_C)+PFid%Cl(M%lmin+i-1,C_C))*M%T(:,M%lmin+i-1)
Cov%C(nl+i,1:nl) = Cov%C(1:nl,2*nl+i)
!TE
Cov%C(1:nl,2*nl+i) = 2*PFid%Cl(M%lmin+i-1,C_C)*M%T(:,M%lmin+i-1)*PFid%Cl(M%lmin:M%lmax,C_C)
Cov%C(2*nl+i,1:nl) = Cov%C(1:nl,2*nl+i)
!XE
Cov%C(nl+1:2*nl,2*nl+i) = rootT(M%lmin+i-1)*rootT(M%lmin:M%lmax)* &
(PFid%Cl(M%lmin:M%lmax,C_C)+PFid%Cl(M%lmin+i-1,C_C))*M%T(:,M%lmin+i-1)
Cov%C(2*nl+i,nl+1:2*nl) = Cov%C(2*nl+i,nl+1:2*nl)
if (vec_size>3) then
!BB
Cov%C(3*nl+1:4*nl,3*nl+i) = 2*PFid%Cl(M%lmin+i-1,C_E)*M%EB(:,M%lmin+i-1)*PFid%Cl(M%lmin:M%lmax,C_E)+&
2 *PFid%Cl(M%lmin+i-1,C_B)*M%EE(:,M%lmin+i-1)*PFid%Cl(M%lmin:M%lmax,C_B)
!EB
Cov%C(3*nl+1:4*nl,2*nl+i) =(rootE(M%lmin+i-1)*rootE(M%lmin:M%lmax)&
+rootB(M%lmin+i-1)*rootB(M%lmin:M%lmax))**2/2*M%EB(:,M%lmin+i-1)
Cov%C(2*nl+i,3*nl+1:4*nl) = Cov%C(3*nl+1:4*nl,2*nl+i)
end if
end do
Minv(nl+1:2*nl,nl+1:2*nl) = M%invX
Minv(2*nl+1:3*nl,2*nl+1:3*nl) = M%invEE
if (vec_size>3) then
Minv(3*nl+1:4*nl,3*nl+1:4*nl) = M%invEE
Minv(2*nl+1:3*nl,3*nl+1:4*nl) = M%invEB
Minv(3*nl+1:4*nl,2*nl+1:3*nl) = M%invEB
end if
else
if (vec_size/=1) stop 'PseudoCl_GetCovariance: vec_size must be 1 for no pol'
end if
print *,'doing mult'
allocate(tmp(nl*vec_size,nl*vec_size))
call Matrix_Mult(Minv,Cov%C,tmp)
call Matrix_Mult_NT(tmp,Minv,Cov%C)
deallocate(tmp)
deallocate(Minv)
print *,'done cov'
end subroutine PseudoCl_GetSimpleCovariance
function sym_ix(n,ax,ay)
integer sym_ix
integer, intent(in) :: n,ax,ay
integer ix, x,y
integer x2,y2
!there's easy analytical result for this too
if (ay>ax) then
x=ay
y=ax
else
x=ax
y=ay
end if
ix=0
do x2=1,n
do y2=1,x2
ix=ix+1
if (x2==x .and. y2==y) then
sym_ix =ix
return
end if
end do
end do
write (*,*) 'n,x,y = ', n,x,y
call MpiStop('sym_ix: error')
end function sym_ix
subroutine TCovMatPolArray_Free(CArr)
Type(TCovMatPolArray) :: CArr
integer i,j,polx,poly, status
do polx=1,CArr%vec_size
do poly =1, polx
call TCovMatArray_Free(CArr%Pol(polx,poly))
end do
end do
deallocate(CArr%Pol)
end subroutine TCovMatPolArray_Free
subroutine EffectiveCovPower(P)
Type(HealpixPower) :: P
real(dp), allocatable :: trueP(:)
real(dp) :: win, norm
integer l, l2, avrange,lmin,lmax
print *,'using effective smoothed C_l'
lmin = 0
lmax= P%lmax
! get effective average C_l for covariance; try just using first one where cosmic var most important
! this tries to cure underestimation of cov because C_l falling rapidly, esp by first trough
!call HealpixPower_Write(SmoothedP(channel),'outfiles/tmpin'//trim(INtToStr(channel)))
avrange=40
allocate(trueP(lmin:lmax))
trueP(lmin:Lmax) = P%Cl(lmin:lmax,C_T)
do l= max(200,avrange+lmin), lmax-avrange
norm = 0
P%Cl(l,C_T)=0
do l2=l-avrange,min(lmax,l+avrange)
win = exp(-real(l2-l)**2/(2*25**2)) !l2^2 and P^2 is OKish
norm = norm + win
P%Cl(l,C_T)= P%Cl(l,C_T) + trueP(l2)*win
end do
P%Cl(l,C_T)=(P%Cl(l,C_T)/norm)
end do
! call HealpixPower_Write(P, 'outfiles/tmp')
deallocate(trueP)
end subroutine EffectiveCovPower
function sym_ix_check(chk, n,i,j) result(res)
logical chk(:)
integer, intent(in) :: n,i,j
integer res
res = sym_ix(n,i,j)
chk(res) = .true.
end function sym_ix_check
subroutine PseudoCl_GetFullCovariance(Coupler, Xi, CArr, PFid, vec_size, Channels, nweights, &
incnoise,incsignal, frac_pt_src_error, pol_weights, step, nstep)
!General result
!At l>30 seems to be very accurate to use just Xi%T rather than Xi%E (and consistent at this approx anyway)
use MatrixUtils
Type(TCovMatPolArray), intent(inout), target :: CArr
logical, intent(in) :: incnoise, incsignal
integer, intent(in) :: step, nstep !to split up so don't blow memory
real, intent(in) :: frac_pt_src_error
logical, intent(in) :: pol_weights
Type(TCouplingMatrix), intent(in) :: Xi(:) !2*ncl(2*ncl+1)/2 vector
Type(TCouplingMatrix), intent(in) :: Coupler(:)
Type(HealpixPower), intent(in) :: PFid
type(TChannel) :: Channels(:)
Type(TCovMat), pointer :: C, EE, EB, BE, BB
integer, intent(in) :: vec_size
integer, intent(in) :: nweights
real(dp) :: ENoiseFac_11,ENoiseFac_22, ENoiseFac_12,ENoiseFac_21,ENoiseFac2
integer aix,nl, ncl_tot, ncl, l, j, status
real fac, crossNoiseFac
real, dimension(:,:), allocatable :: rootT, rootE, rootB
real, dimension(:), allocatable :: rootT_11,RootT_12,rootT_21,rootT_22
real, dimension(:), allocatable :: BeamAvTE, rootE_11,RootE_12,rootE_21,rootE_22
real, dimension(:), allocatable :: rootB_11,RootB_12,rootB_21,rootB_22
real(dp), dimension(:,:), allocatable :: tmp, tmp11,tmp12,tmp21,tmp22
integer, parameter :: w2=1, wTT=2, w2T =3
integer i1,i2,j1,j2,ix, x,y, ix2, x2,y2, lmin, lmax
integer ix_couple,ix2_couple
integer s_11_s_22,s_12_s_21, n_11_s_22,n_12_s_21, s_11_n_22,s_12_n_21, n_11_n_22,n_12_n_21
integer PP_s_11_s_22, PP_s_12_s_21, PP_n_11_s_22, PP_n_12_s_21, PP_s_11_n_22, PP_s_12_n_21, PP_n_11_n_22, PP_n_12_n_21
integer PT_s_11_s_22, PT_s_12_s_21, PT_n_11_s_22, PT_n_12_s_21, PT_s_11_n_22, PT_s_12_n_21, PT_n_11_n_22, PT_n_12_n_21
integer MM_s_11_s_22, MM_s_12_s_21
integer MT_s_11_n_22, MT_s_12_n_21, MP_s_12_n_21, PM_n_11_s_22
integer MT_s_11_s_22, MT_s_12_s_21, PM_s_11_s_22, PM_s_12_s_21
logical dopol
integer polx,poly, nchannels, channel
double precision STime
Type(HealpixPower), allocatable :: DummyNoise(:)
Type(HealpixPower), allocatable :: SmoothedP(:)
Type(HealpixPower) :: AvP
Type(TBeam) :: AvBeam
integer chanx, chanx2,chany,chany2
integer nvarmaps, nmaskmaps
integer ix_11,ix_12, ix_21,ix_22
logical, allocatable :: chk(:)
nchannels = size(Channels)
lmin = Coupler(1)%lmin
lmax = Coupler(1)%lmax
ncl = (nweights*(nweights+1))/2
ncl_tot = (nweights*nchannels*(nweights*nchannels+1))/2
print *,'PseudoCl_GetFullCovariance: channels = ',nchannels
dopol= Coupler(1)%has_pol
nl = lmax-lmin+1
if (pol_weights) then
nvarmaps= (2*nchannels+3)*ncl
nmaskmaps=(nchannels+1)*ncl
else
nvarmaps= (nchannels+1)*ncl
nmaskmaps =0
end if
allocate(chk(size(Xi)))
chk = .false.
if (step /= 0) then
allocate(rootT_11(lmin:lmax),rootT_21(lmin:lmax),rootT_12(lmin:lmax),rootT_22(lmin:lmax))
allocate(rootT(lmin:lmax,nchannels))
if (dopol) then
allocate(rootE_11(lmin:lmax),rootE_21(lmin:lmax),rootE_12(lmin:lmax),rootE_22(lmin:lmax))
allocate(rootB_11(lmin:lmax),rootB_21(lmin:lmax),rootB_12(lmin:lmax),rootB_22(lmin:lmax))
allocate(RootE(lmin:lmax,nchannels),RootB(lmin:lmax,nchannels))
allocate(BeamAvTE(lmin:lmax))
end if
call HealpixPower_Nullify(AvP)
print*,'Noise Scale= ', NoiseScale
crossNoiseFac = NoiseScale**2*cross_noise_scale ! for the N-N terms, the scaling -can be increased for cross-spectra
if (.not. associated(Coupler(1)%InvT)) call MpiStop('Need inverses')
print *,'getting smoothed Cls'
allocate(SmoothedP(nchannels))
do channel = 1, nchannels
call healpixPower_Nullify(SmoothedP(channel))
call HealpixPower_Assign(SmoothedP(channel), PFid)
if (incsignal) then
SmoothedP(channel)%Cl(:,C_T) = SmoothedP(channel)%Cl(:,C_T) + Channels(channel)%PtSrcA
call TBeam_PowerSmooth(Channels(channel)%Beam, SmoothedP(channel), -1)
call EffectiveCovPower(SmoothedP(channel))
else
SmoothedP(channel)%Cl=0
end if
rootT(:,channel) = sqrt(SmoothedP(channel)%Cl(lmin:lmax,C_T))
if (dopol) then
! if (nchannels >1) call MpiStop('not done multi-channel pol')
rootE(:,channel) = sqrt(SmoothedP(channel)%Cl(lmin:lmax,C_E))