-
Notifications
You must be signed in to change notification settings - Fork 7
/
cmbmain.f90
2580 lines (2088 loc) · 89.9 KB
/
cmbmain.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
! This this is the main CAMB program module.
!
! Code for Anisotropies in the Microwave Background
! by Antony lewis (http://cosmologist.info) and Anthony Challinor
! See readme.html for documentation.
! Note that though the code is internally parallelised, it is not thread-safe
! so you cannot generate more than one model at the same time in different threads.
!
! Based on CMBFAST by Uros Seljak and Matias Zaldarriaga, itself based
! on Boltzmann code written by Edmund Bertschinger, Chung-Pei Ma and Paul Bode.
! Original CMBFAST copyright and disclaimer:
!
! Copyright 1996 by Harvard-Smithsonian Center for Astrophysics and
! the Massachusetts Institute of Technology. All rights reserved.
!
! THIS SOFTWARE IS PROVIDED "AS IS", AND M.I.T. OR C.f.A. MAKE NO
! REPRESENTATIONS OR WARRANTIES, EXPRESS OR IMPlIED.
! By way of example, but not limitation,
! M.I.T. AND C.f.A MAKE NO REPRESENTATIONS OR WARRANTIES OF
! MERCHANTABIlITY OR FITNESS FOR ANY PARTICUlAR PURPOSE OR THAT
! THE USE OF THE lICENSED SOFTWARE OR DOCUMENTATION WIll NOT INFRINGE
! ANY THIRD PARTY PATENTS, COPYRIGHTS, TRADEMARKS OR OTHER RIGHTS.
!
! portions of this software are based on the COSMICS package of
! E. Bertschinger. See the lICENSE file of the COSMICS distribution
! for restrictions on the modification and distribution of this software.
module CAMBmain
! This code evolves the linearized perturbation equations of general relativity,
! the Boltzmann equations and the fluid equations for perturbations
! of a Friedmann-Robertson-Walker universe with a supplied system of gauge-dependent equation
! in a modules called GaugeInterface. The sources for the line of sight integral are
! computed at sampled times during the evolution for various of wavenumbers. The sources
! are then interpolated to a denser wavenumber sampling for computing the line of
! sight integrals of the form Integral d(conformal time) source_k * bessel_k_l.
! For CP%flat models the bessel functions are interpolated from a pre-computed table, for
! non-CP%flat models the hyperspherical Bessel functions are computed by integrating their
! differential equation. Both phases ('Evolution' and 'Integration') can do separate
! wavenumbers in parallel.
! The time variable is conformal time dtau=dt/a(t) and the spatial dependence is Fourier transformed
! with q=sqrt(k**2 + (|m|+1)K), comoving distances are x=CP%r/a(t), with a(t)=1 today.
! The units of both length and time are Mpc.
! Many elements are part of derived types (to make thread safe or to allow non-sequential code use
! CP = CAMB parameters
! EV = Time evolution variables
! IV = Source integration variables
! CT = Cl transfer data
! MT = matter transfer data
! Modules are defined in modules.f90, except GaugeInterface which gives the background and gauge-dependent
! perturbation equations, and InitialPower which provides the initial power spectrum.
use precision
use ModelParams
use ModelData
use GaugeInterface
use Transfer
use SpherBessels
use lvalues
use MassiveNu
use InitialPower
use Errors
implicit none
private
logical :: WantLateTime = .false. !if lensing or redshift windows
logical ExactClosedSum !do all nu values in sum for Cls for Omega_k>0.1
!Variables for integrating the sources with the bessel functions for each wavenumber
type IntegrationVars
integer q_ix
real(dl) q, dq !q value we are doing and delta q
! real(dl), dimension(:,:), pointer :: Delta_l_q
!Contribution to C_l integral from this k
real(dl), dimension(:,:), pointer :: Source_q, ddSource_q
!Interpolated sources for this k
integer SourceSteps !number of steps up to where source is zero
end type IntegrationVars
integer SourceNum
!SourceNum is total number sources (2 or 3 for scalars, 3 for tensors).
real(dl) tautf(0:max_transfer_redshifts) !Time of Trasfer%redshifts
real(dl), dimension(:,:,:), allocatable :: Src, ddSrc !Sources and second derivs
! indices Src( k_index, source_index, time_step_index )
real(dl), dimension(:,:,:), allocatable :: iCl_scalar, iCl_vector, iCl_tensor
! Cls at the l values we actually compute, iCl_xxx(l_index, Cl_type, initial_power_index)
real(dl), dimension(:,:,:,:), allocatable :: iCl_Array
!Full covariance at each L (alternative more general arrangement to above)
! values of q to evolve the propagation equations to compute the sources
Type(Regions) :: Evolve_q
real(dl),parameter :: qmin0=0.1_dl
real(dl) :: dtaurec_q
! qmax - CP%Max_eta_k/CP%tau0, qmin = qmin0/CP%tau0 for flat case
real(dl) qmin, qmax
real(dl) max_etak_tensor , max_etak_vector, max_etak_scalar
! Will only be calculated if k*tau < max_etak_xx
integer maximum_l !Max value of l to compute
real(dl) :: maximum_qeta = 3000._dl
integer :: l_smooth_sample = 3000 !assume transfer functions effectively small for k>2*l_smooth_sample
real(dl) :: fixq = 0._dl !Debug output of one q
real(dl) :: ALens = 1._dl
Type(ClTransferData), pointer :: ThisCT
public cmbmain, ALens, ClTransferToCl, InitVars !InitVars for BAO hack
contains
subroutine cmbmain
integer q_ix
! real clock_start, clock_stop ! RH timing
! real clock_totstart, clock_totstop ! RH timing
type(EvolutionVars) EV
! Timing variables for testing purposes. Used if DebugMsgs=.true. in ModelParams
real(sp) actual,timeprev,starttime
! clock_totstart = 0
! call cpu_time(clock_totstart) ! RH timing
WantLateTime = CP%DoLensing .or. num_redshiftwindows > 0
if (CP%WantCls) then
if (CP%WantTensors .and. CP%WantScalars) stop 'CMBMAIN cannot generate tensors and scalars'
!Use CAMB_GetResults instead
if (CP%WantTensors) then
maximum_l = CP%Max_l_tensor
maximum_qeta = CP%Max_eta_k_tensor
else
maximum_l = CP%Max_l
maximum_qeta = CP%Max_eta_k
end if
! call cpu_time(clock_start) ! RH timing
call initlval(lSamp, maximum_l)
! call cpu_time(clock_stop) ! RH timing
! print*, 'this is how long initlval takes in cmbmain', clock_stop - clock_start
end if
if (DebugMsgs .and. Feedbacklevel > 0) then
actual=GetTestTime()
starttime=actual !times don't include reading the Bessel file
end if
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call InitVars !Most of single thread time spent here (in InitRECFAST)
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'this is how long it takes to InitVars in cmbmain', clock_stop - clock_start
if (global_error_flag/=0) return
if (DebugMsgs .and. Feedbacklevel > 0) then
timeprev=actual
actual=GetTestTime()
write(*,*) actual-timeprev,' Timing for InitVars'
write (*,*) 'r = ',real(CP%r),' scale = ',real(scale), 'age = ', real(CP%tau0)
end if
!JD 08/13 for nonlinear lensing of CMB + MPK compatibility
!if (.not. CP%OnlyTransfers .or. CP%NonLinear==NonLinear_Lens) call InitializePowers(CP%InitPower,CP%curv)
if (.not. CP%OnlyTransfers .or. CP%NonLinear==NonLinear_Lens .or. CP%NonLinear==NonLinear_both) then
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call InitializePowers(CP%InitPower,CP%curv)
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing of InitializePowers in cmbmain', clock_stop - clock_start
end if
if (global_error_flag/=0) return
!Calculation of the CMB sources.
if (CP%WantCls) then
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call SetkValuesForSources
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for SetkValuesForSources in cmbmain', clock_stop - clock_start
end if
if (CP%WantTransfer) then
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call InitTransfer
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for InitTransfer in cmbmain', clock_stop - clock_start
end if
!***note that !$ is the prefix for conditional multi-processor compilation***
!$ if (ThreadNum /=0) call OMP_SET_NUM_THREADS(ThreadNum)
if (CP%WantCls) then
if (DebugMsgs .and. Feedbacklevel > 0) write(*,*) 'Set ',Evolve_q%npoints,' source k values'
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call GetSourceMem
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for GetSourceMem in cmbmain', clock_stop - clock_start
if (CP%WantScalars) then
ThisCT => CTransScal
else if (CP%WantVectors) then
ThisCT => CTransVec
else
ThisCT => CTransTens
end if
ThisCT%NumSources = SourceNum
ThisCT%ls = lSamp
!$OMP PARAllEl DO DEFAUlT(SHARED),SCHEDUlE(DYNAMIC) &
!$OMP & PRIVATE(EV, q_ix)
! print*, 'Im going to DoSourcek for Evolve_q%npoints times: ', Evolve_q%npoints
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
do q_ix= 1,Evolve_q%npoints
if (global_error_flag==0) then
call DoSourcek(EV,q_ix)
end if
end do
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for DoSourceK in cmbmain', clock_stop - clock_start
!$OMP END PARAllEl DO
if (DebugMsgs .and. Feedbacklevel > 0) then
timeprev=actual
actual=GetTestTime()
write(*,*) actual-timeprev,' Timing for source calculation'
end if
endif !WantCls
! If transfer functions are requested, set remaining k values and output
if (CP%WantTransfer .and. global_error_flag==0) then
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call TransferOut
if (DebugMsgs .and. Feedbacklevel > 0) then
timeprev=actual
actual=GetTestTime()
write(*,*) actual-timeprev,' Timing for transfer k values'
end if
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for TransferOut in cmbmain', clock_stop - clock_start
end if
if (CP%WantTransfer .and. CP%WantCls .and. WantLateTime &
.and. (CP%NonLinear==NonLinear_Lens .or. CP%NonLinear==NonLinear_both) .and. global_error_flag==0) then
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call MakeNonlinearSources
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for MakeNonlinearSources in cmbmain', clock_stop - clock_start
if (DebugMsgs .and. Feedbacklevel > 0) then
timeprev=actual
actual=GetTestTime()
write(*,*) actual-timeprev,' Timing for NonLinear sources'
end if
end if
if (CP%WantTransfer .and. .not. CP%OnlyTransfers .and. global_error_flag==0) then
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call Transfer_Get_sigma8(MT,8._dl)
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for Transfer_Get_sigma8 in cmbmain', clock_stop - clock_start
end if
!Can call with other arguments if need different size
! if CMB calculations are requested, calculate the Cl by
! integrating the sources over time and over k.
if (CP%WantCls) then
if (global_error_flag==0) then
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call InitSourceInterpolation
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for InitSourceInterpolation in cmbmain', clock_stop - clock_start !RH timing
ExactClosedSum = CP%curv > 5e-9_dl .or. scale < 0.93_dl
max_bessels_l_index = ThisCT%ls%l0
max_bessels_etak = maximum_qeta
if (CP%WantScalars) then
! clock_start =0
! call cpu_time(clock_start) ! RH timing
call GetLimberTransfers
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for GetLimberTransfers in cmbmain', clock_stop - clock_start
end if
ThisCT%max_index_nonlimber = max_bessels_l_index
if (CP%flat) call InitSpherBessels
!This is only slow if not called before with same (or higher) Max_l, Max_eta_k
!Preferably stick to Max_l being a multiple of 50
call SetkValuesForInt
if (DebugMsgs .and. Feedbacklevel > 0) write(*,*) 'Set ',ThisCT%q%npoints,' integration k values'
!Begin k-loop and integrate Sources*Bessels over time
!$OMP PARAllEl DO DEFAUlT(SHARED),SHARED(TimeSteps), SCHEDUlE(STATIC,4)
! print*, 'SourceToTransfers is going to integrate ThisCT%q%npoints times', ThisCT%q%npoints
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
do q_ix=1,ThisCT%q%npoints
call SourceToTransfers(q_ix)
end do !q loop
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for SourceToTransfersLoop in cmbmain', clock_stop - clock_start ! RH timing
!$OMP END PARAllEl DO
if (DebugMsgs .and. Feedbacklevel > 0) then
timeprev=actual
actual=GetTestTime()
write(*,*)actual-timeprev,' Timing For Integration'
end if
end if
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
call FreeSourceMem
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'timing for FreeSourceMem in cmbmain', clock_stop - clock_start
!Final calculations for CMB output unless want the Cl transfer functions only.
if (.not. CP%OnlyTransfers .and. global_error_flag==0) then
! clock_start = 0
! call cpu_time(clock_start) ! RH timing
!!! RH axions need to be called again here !!!
call ClTransferToCl(CTransScal,CTransTens, CTransVec)
! clock_stop = 0
! call cpu_time(clock_stop) ! RH timing
! print*, 'this is the timing for ClTransferToCl in cmbmain', clock_stop - clock_start
end if
end if
if (DebugMsgs .and. Feedbacklevel > 0) then
timeprev=actual
actual = GetTestTime()
write(*,*) actual - timeprev,' Timing for final output'
write(*,*) actual -starttime,' Timing for whole of cmbmain'
end if
! clock_totstop = 0
! call cpu_time(clock_totstop) ! RH timing
! print*, 'timing for the whole of cmbmain inside it', clock_totstop - clock_totstart ! RH timing
end subroutine cmbmain
subroutine ClTransferToCl(CTransS,CTransT, CTransV)
Type(ClTransferData) :: CTransS,CTransT, CTransV
if (CP%WantScalars .and. global_error_flag==0) then
lSamp = CTransS%ls
allocate(iCl_Scalar(CTransS%ls%l0,C_Temp:C_last,CP%InitPower%nn))
iCl_scalar = 0
if (has_cl_2D_array) then
allocate(iCl_Array(CTransS%ls%l0,CTransS%NumSources,CTransS%NumSources,CP%InitPower%nn))
iCl_Array = 0
end if
call CalcLimberScalCls(CTransS)
call CalcScalCls(CTransS)
if (DebugMsgs .and. Feedbacklevel > 0) write (*,*) 'CalcScalCls'
end if
if (CP%WantVectors .and. global_error_flag==0) then
allocate(iCl_vector(CTransV%ls%l0,C_Temp:CT_Cross,CP%InitPower%nn))
iCl_vector = 0
call CalcVecCls(CTransV,GetInitPowerArrayVec)
if (DebugMsgs .and. Feedbacklevel > 0) write (*,*) 'CalcVecCls'
end if
if (CP%WantTensors .and. global_error_flag==0) then
allocate(iCl_Tensor(CTransT%ls%l0,CT_Temp:CT_Cross,CP%InitPower%nn))
iCl_tensor = 0
call CalcTensCls(CTransT,GetInitPowerArrayTens)
if (DebugMsgs .and. Feedbacklevel > 0) write (*,*) 'CalcTensCls'
end if
if (global_error_flag==0) then
call Init_Cls
! Calculating Cls for every l.
call InterpolateCls(CTransS,CTransT, CTransV)
if (DebugMsgs .and. Feedbacklevel > 0) write (*,*) 'InterplolateCls'
end if
if (CP%WantScalars .and. allocated(iCl_Scalar)) deallocate(iCl_scalar)
if (CP%WantScalars .and. allocated(iCl_Array)) deallocate(iCl_Array)
if (CP%WantVectors .and. allocated(iCl_Vector)) deallocate(iCl_vector)
if (CP%WantTensors .and. allocated(iCl_Tensor)) deallocate(iCl_tensor)
if (global_error_flag/=0) return
if (CP%OutputNormalization >=2) call NormalizeClsAtl(CP%OutputNormalization)
!Normalize to C_l=1 at l=OutputNormalization
end subroutine ClTransferToCl
subroutine CalcLimberScalCls(CTrans)
Type(ClTransferData) :: CTrans
integer ell, i, s_ix, pix
real(dl) CL, reall,fac, dbletmp
integer s_ix2,j,n
Type(LimberRec), pointer :: LimbRec,LimbRec2
integer winmin
if (limber_phiphi>0) then
winmin = 0
else
winmin=1
end if
do pix=1,CP%InitPower%nn
do i =winmin, num_redshiftwindows
s_ix = 3+i
if (CTrans%limber_l_min(s_ix) /=0) then
do j= i, num_redshiftwindows
s_ix2 = 3+j
if (CTrans%limber_l_min(s_ix2) /=0) then
!$OMP PARALLEL DO DEFAUlT(SHARED), SCHEDUlE(STATIC,2), PRIVATE(Cl,ell,reall,fac,dbletmp,n,LimbRec,LimbRec2)
do ell = max(CTrans%limber_l_min(s_ix), CTrans%limber_l_min(s_ix2)), Ctrans%ls%l0
Cl = 0
reall = real(CTrans%ls%l(ell),dl)
fac = (2*pi**2)/fourpi/(reall+0.5_dl)**3 !fourpi because multipled by fourpi later
LimbRec => CTrans%Limber_windows(s_ix,ell)
LimbRec2 => CTrans%Limber_windows(s_ix2,ell)
do n = max(LimbRec%n1,LimbRec2%n1), min(LimbRec%n2,LimbRec2%n2)
!Actually integral over chi; source has sqrt( chi dchi)
!Same n corresponds to same k since ell fixed here
Cl = Cl + LimbRec%Source(n)*LimbRec2%Source(n) * fac * ScalarPower(LimbRec%k(n) ,pix)
end do
if(j==0 .and. i==0) iCl_scalar(ell,C_Phi,pix) = Cl
if (has_cl_2D_array) then
dbletmp=(reall*(reall+1))/OutputDenominator*fourpi
iCl_Array(ell,s_ix,s_ix2,pix) = Cl*dbletmp
if (i/=j) iCl_Array(ell,s_ix2,s_ix,pix)=iCl_Array(ell,s_ix,s_ix2,pix)
end if
end do
!$OMP END PARAllEl DO
end if
end do
end if
end do
end do
end subroutine CalcLimberScalCls
subroutine GetLimberTransfers
integer ell, ell_needed
integer i,s_ix
integer n1,n2,n, ell_limb
real(dl) int,k, chi, chimin
integer klo, khi
real(dl) a0,b0,ho2o6,a03,b03,ho
Type(LimberRec), pointer :: LimbRec
call Init_Limber(ThisCT)
if (num_redshiftwindows==0 .and. limber_phiphi==0) return
if (ThisCT%ls%l(ThisCT%ls%l0) > 5000) then
max_bessels_l_index = lvalues_indexOf(ThisCT%ls,5000)
else
max_bessels_l_index = ThisCT%ls%l0
end if
if (CP%Want_CMB) then
max_bessels_etak= min(ThisCT%ls%l(ThisCT%ls%l0),3000)*2.5_dl*AccuracyBoost
else
max_bessels_etak = 5000
end if
do i =0, num_redshiftwindows
s_ix = 3+i
if (i==0) then
ell_limb = limber_phiphi
end if
ell_needed = ThisCT%ls%l(ThisCT%ls%l0)
do ell = 1, ThisCT%ls%l0
if (ThisCT%ls%l(ell) >= ell_limb) then
ThisCT%limber_l_min(s_ix) = ell
ell_needed = ThisCT%ls%l(ell)
max_bessels_l_index = max(max_bessels_l_index,ThisCT%limber_l_min(s_ix)-1)
if (FeedbackLevel > 1) write (*,*) i,'Limber switch', ell_needed
exit
end if
end do
if (ThisCT%limber_l_min(s_ix)/=0) then
if (i==0) then
n1 = Ranges_IndexOf(TimeSteps,tau_maxvis)
n2 = TimeSteps%npoints-1
end if
do ell = ThisCT%limber_l_min(s_ix), ThisCT%ls%l0
LimbRec => ThisCT%Limber_windows(s_ix,ell)
LimbRec%n1 = n1
LimbRec%n2 = n2
allocate(LimbRec%k(n1:n2))
allocate(LimbRec%Source(n1:n2))
int = 0
do n = n1,n2
chi = (CP%tau0-TimeSteps%points(n))
k = (ThisCT%ls%l(ell)+0.5_dl)/chi
LimbRec%k(n) = k
if (k<=qmax) then
klo = Ranges_IndexOf(Evolve_q, k)
khi=klo+1
ho=Evolve_q%points(khi)-Evolve_q%points(klo)
a0=(Evolve_q%points(khi)-k)/ho
b0=(k-Evolve_q%points(klo))/ho
ho2o6 = ho**2/6
a03=(a0**3-a0)
b03=(b0**3-b0)
LimbRec%Source(n)= sqrt(chi*TimeSteps%dpoints(n))* (a0*Src(klo,s_ix,n)+&
b0*Src(khi,s_ix,n)+(a03 *ddSrc(klo,s_ix,n)+ b03*ddSrc(khi,s_ix,n)) *ho2o6)
else
LimbRec%Source(n)=0
end if
end do
end do
else
max_bessels_l_index = ThisCT%ls%l0
end if
end do
end subroutine GetLimberTransfers
subroutine SourceToTransfers(q_ix)
integer q_ix
type(IntegrationVars) :: IV
allocate(IV%Source_q(TimeSteps%npoints,SourceNum))
if (.not.CP%flat) allocate(IV%ddSource_q(TimeSteps%npoints,SourceNum))
call IntegrationVars_init(IV)
IV%q_ix = q_ix
IV%q =ThisCT%q%points(q_ix)
IV%dq= ThisCT%q%dpoints(q_ix)
call InterpolateSources(IV)
call DoSourceIntegration(IV)
if (.not.CP%flat) deallocate(IV%ddSource_q)
deallocate(IV%Source_q)
end subroutine SourceToTransfers
subroutine InitTransfer
integer nu,lastnu, ntodo, nq, q_ix, first_i
real(dl) dlog_lowk1,dlog_lowk, d_osc,dlog_osc, dlog_highk, boost
real(dl) amin,q_switch_lowk,q_switch_lowk1,q_switch_osc,q_switch_highk
real(dl), dimension(:), allocatable :: q_transfer
! print*, 'In InitTransfer the kmax value being used is: ', CP%Transfer%kmax
if (CP%Transfer%k_per_logint==0) then
!Optimized spacing
!Large log spacing on superhorizon scales
!Linear spacing for horizon scales and first few baryon osciallations
!Log spacing for last few oscillations
!large log spacing for small scales
boost = AccuracyBoost
if (CP%Transfer%high_precision) boost = boost*1.5
q_switch_lowk1 = 0.7/taurst
dlog_lowk1=2*boost
q_switch_lowk = 8/taurst
dlog_lowk=8*boost
if (HighAccuracyDefault) dlog_lowk = dlog_lowk*2.5
q_switch_osc = min(CP%Transfer%kmax,30/taurst)
d_osc= 200*boost
if (HighAccuracyDefault) d_osc = d_osc*1.8
q_switch_highk = min(CP%Transfer%kmax,60/taurst)
dlog_osc = 17*boost
if (HighAccuracyDefault) q_switch_highk = min(CP%Transfer%kmax,90/taurst)
!Then up to kmax
dlog_highk = 3*boost
amin = 5e-5_dl
nq=int((log(CP%Transfer%kmax/amin))*d_osc)+1
allocate(q_transfer(nq))
nq=int((log(q_switch_lowk1/amin))*dlog_lowk1)+1
do q_ix=1, nq
q_transfer(q_ix) = amin*exp((q_ix-1)/dlog_lowk1)
end do
MT%num_q_trans = nq
nq=int(log( q_switch_lowk/q_transfer(MT%num_q_trans))*dlog_lowk) +1
do q_ix=1, nq
q_transfer(MT%num_q_trans+q_ix) = q_transfer(MT%num_q_trans)*exp(q_ix/dlog_lowk)
end do
MT%num_q_trans = MT%num_q_trans + nq
nq=int((q_switch_osc-q_transfer(MT%num_q_trans))*d_osc)+1
do q_ix=1, nq
q_transfer(MT%num_q_trans+q_ix) = q_transfer(MT%num_q_trans)+ q_ix/d_osc
end do
MT%num_q_trans = MT%num_q_trans + nq
if (CP%Transfer%kmax > q_transfer(MT%num_q_trans)) then
nq=int(log( q_switch_highk/q_transfer(MT%num_q_trans))*dlog_osc) +1
do q_ix=1, nq
q_transfer(MT%num_q_trans+q_ix) = q_transfer(MT%num_q_trans)*exp(q_ix/dlog_osc)
end do
MT%num_q_trans = MT%num_q_trans + nq
end if
if (CP%Transfer%kmax > q_transfer(MT%num_q_trans)) then
nq=int(log(CP%Transfer%kmax/q_transfer(MT%num_q_trans))*dlog_highk)+1
do q_ix=1, nq
q_transfer(MT%num_q_trans+q_ix) = q_transfer(MT%num_q_trans)*exp(q_ix/dlog_highk)
end do
MT%num_q_trans = MT%num_q_trans + nq
end if
else
!Fixed spacing
! MT%num_q_trans=int((log(CP%Transfer%kmax)-log(qmin))*CP%Transfer%k_per_logint)+1
! allocate(q_transfer(MT%num_q_trans))
! do q_ix=1, MT%num_q_trans
! q_transfer(q_ix) = qmin*exp(real(q_ix)/CP%Transfer%k_per_logint)
! end do
MT%num_q_trans=int((log(CP%Transfer%kmax)-log(qmin))*CP%Transfer%k_per_logint)+1
allocate(q_transfer(MT%num_q_trans))
do q_ix=1, MT%num_q_trans
q_transfer(q_ix) = qmin*exp(real(q_ix)/CP%Transfer%k_per_logint)
end do
end if
if (CP%closed) then
lastnu=0
ntodo = 0
do q_ix=1,MT%num_q_trans
nu =nint(CP%r*q_transfer(q_ix))
if (.not. ((nu<3).or.(nu<=lastnu))) then
ntodo=ntodo+1
q_transfer(ntodo)= nu/CP%r
lastnu=nu
end if
end do
MT%num_q_trans = ntodo
end if
if (CP%WantCls) then
ntodo = MT%num_q_trans
first_i = ntodo+1
do q_ix = 1,ntodo
if (q_transfer(q_ix) > Evolve_q%highest+1d-4) then
!Feb13 fix for case with closed universe where qmax is not neccessarily right quantized value
! if (q_transfer(q_ix) > qmax) then
first_i=q_ix
exit
end if
end do
if (first_i > ntodo) then
MT%num_q_trans = Evolve_q%npoints
else
MT%num_q_trans = Evolve_q%npoints + (ntodo - first_i+1)
end if
call Transfer_Allocate(MT)
MT%q_trans(1:Evolve_q%npoints) = Evolve_q%points(1:Evolve_q%npoints)
if (MT%num_q_trans > Evolve_q%npoints) then
MT%q_trans(Evolve_q%npoints+1:MT%num_q_trans) = q_transfer(first_i:ntodo)
end if
else
Evolve_q%npoints = 0
call Transfer_Allocate(MT)
MT%q_trans = q_transfer(1:MT%num_q_trans)
end if
deallocate(q_transfer)
end subroutine InitTransfer
function GetTauStart(q)
real(dl), intent(IN) :: q
real(dl) taustart, GetTauStart
real(dl) tauosc,om,rhomass,grhonu,taurad,tauax_init,taueq
double precision arad,abar,taubar
! Begin when wave is far outside horizon.
! DM: also begin when tau<tau_osc
! Conformal time (in Mpc) in the radiation era, for photons plus 3 species
! of relativistic neutrinos.
if (CP%flat) then
taustart=0.001_dl/q
else
taustart=0.001_dl/sqrt(q**2-CP%curv)
end if
if (fixq/=0._dl) then
taustart = 0.001_dl/fixq
end if
! Make sure to start early in the radiation era.
taustart=min(taustart,0.1_dl)
! Start when massive neutrinos are strongly relativistic.
if (CP%Num_nu_massive>0) then
taustart=min(taustart,1.d-3/maxval(nu_masses(1:CP%Nu_mass_eigenstates))/adotrad)
end if
! Start when axions are not oscillating
!update start time to use a_init from w_evolve routine 9/19 dgrin, comment out other checks (superfluous)
if(CP%Omegaax>0) then
rhomass = sum(grhormass(1:CP%Nu_mass_eigenstates))
grhonu=rhomass+grhornomass
om = (grhob+grhoc)/sqrt(3.0d0*(grhog+grhonu))
tauosc= (a_osc*0.01d0/adotrad)/om
!DM:
arad=((a_osc**3.0d0)*((grhog+grhonu)*0.01d0)/(grhoax))**(1.0d0/4.0d0)
abar=((a_osc**3.0d0)*((grhob)*0.01d0)/(grhoax))**(1.0d0/3.0d0)
taueq= (-1.0d0+dsqrt(1.0d0+CP%aeq*0.01d0/adotrad))/om/2.0d0
!calculate corresponding conformal time
!taurad= (-1.0d0+dsqrt(1.0d0+arad/adotrad))/om/2.0d0
!taubar= (-1.0d0+dsqrt(1.0d0+abar/adotrad))/om/2.0d0
!old working start time for other modes!
taustart=min(taustart,tauosc,taueq)
end if
GetTauStart=taustart
end function GetTauStart
subroutine DoSourcek(EV,q_ix)
integer q_ix
real(dl) taustart
type(EvolutionVars) EV
EV%q=Evolve_q%points(q_ix)
if (fixq/=0._dl) then
EV%q= min(500._dl,fixq) !for testing
end if
EV%q2=EV%q**2
EV%q_ix = q_ix
EV%TransferOnly=.false.
taustart = GetTauStart(EV%q)
if (fixq/=0._dl) then
EV%q= fixq
EV%q2=EV%q**2
end if
call GetNumEqns(EV)
if (CP%WantScalars .and. global_error_flag==0) call CalcScalarSources(EV,taustart)
if (CP%WantVectors .and. global_error_flag==0) call CalcVectorSources(EV,taustart)
if (CP%WantTensors .and. global_error_flag==0) call CalcTensorSources(EV,taustart)
end subroutine DoSourcek
subroutine GetSourceMem
if (CP%WantScalars) then
if (WantLateTime) then
SourceNum=3
C_last = C_PhiE
SourceNum=SourceNum + num_redshiftwindows + num_extra_redshiftwindows
else
SourceNum=2
C_last = C_Cross
end if
else
SourceNum=3
end if
allocate(Src(Evolve_q%npoints,SourceNum,TimeSteps%npoints))
Src=0
allocate(ddSrc(Evolve_q%npoints,SourceNum,TimeSteps%npoints))
end subroutine GetSourceMem
subroutine FreeSourceMem
if (allocated(Src))deallocate(Src, ddSrc)
call Ranges_Free(Evolve_q)
end subroutine FreeSourceMem
! initial variables, number of steps, etc.
subroutine InitVars
use ThermoData
use precision
use ModelParams
implicit none
real(dl) taumin, maxq, initAccuracyBoost
integer itf
initAccuracyBoost = AccuracyBoost
! Maximum and minimum k-values.
if (CP%flat) then
qmax=maximum_qeta/CP%tau0
qmin=qmin0/CP%tau0/initAccuracyBoost
else
qmax=maximum_qeta/CP%r/CP%chi0
qmin=qmin0/CP%r/CP%chi0/initAccuracyBoost
end if
! Timesteps during recombination (tentative, the actual
! timestep is the minimum between this value and taurst/40,
! where taurst is the time when recombination starts - see inithermo
dtaurec_q=4/qmax/initAccuracyBoost
if (.not. CP%flat) dtaurec_q=dtaurec_q/6
!AL:Changed Dec 2003, dtaurec feeds back into the non-flat integration via the step size
dtaurec = dtaurec_q
!dtau rec may be changed by inithermo
max_etak_tensor = initAccuracyBoost*maximum_qeta /10
max_etak_scalar = initAccuracyBoost*max(1700._dl,maximum_qeta) /20
if (maximum_qeta <3500 .and. AccuracyBoost < 2) max_etak_scalar = max_etak_scalar * 1.5
!tweak to get large scales right
max_etak_vector = max_etak_scalar
if (CP%WantCls) then
maxq = qmax
if (CP%WantTransfer) maxq=max(qmax,CP%Transfer%kmax)
else
maxq=CP%Transfer%kmax
end if
taumin=GetTauStart(maxq)
! Initialize baryon temperature and ionization fractions vs. time.
! This subroutine also fixes the timesteps where the sources are
! saved in order to do the integration. So TimeSteps is set here.
!These routines in ThermoData (modules.f90)
call inithermo(taumin,CP%tau0)
if (global_error_flag/=0) return
if (DebugMsgs .and. Feedbacklevel > 0) write (*,*) 'inithermo'
!Do any array initialization for propagation equations
call GaugeInterface_Init
if (Feedbacklevel > 0) &
write(*,'("tau_recomb/Mpc = ",f7.2," tau_now/Mpc = ",f8.1)') tau_maxvis,CP%tau0
! Calculating the times for the outputs of the transfer functions.
!
if (CP%WantTransfer) then
do itf=1,CP%Transfer%num_redshifts
tautf(itf)=min(TimeOfz(CP%Transfer%redshifts(itf)),CP%tau0)
if (itf>1) then
if (tautf(itf) <= tautf(itf-1)) then
stop 'Transfer redshifts not set or out of order'
end if
end if
end do
endif
end subroutine InitVars
subroutine SetkValuesForSources
implicit none
real(dl) dlnk0, dkn1, dkn2, q_switch, q_cmb, dksmooth
real(dl) qmax_log
real(dl) SourceAccuracyBoost
! set k values for which the sources for the anisotropy and
! polarization will be calculated. For low values of k we
! use a logarithmic spacing. closed case dealt with by SetClosedkValues
SourceAccuracyBoost = AccuracyBoost
if (CP%WantScalars .and. CP%Reion%Reionization .and. CP%AccuratePolarization) then
dlnk0=2._dl/10/SourceAccuracyBoost
!Need this to get accurate low l polarization
else
dlnk0=5._dl/10/SourceAccuracyBoost
if (CP%closed) dlnk0=dlnk0/2
end if
if (CP%AccurateReionization) dlnk0 = dlnk0/2
dkn1=0.6_dl/taurst/SourceAccuracyBoost
dkn2=0.9_dl/taurst/SourceAccuracyBoost
if (HighAccuracyDefault) dkn2=dkn2/1.2
if (CP%WantTensors .or. CP%WantVectors) then
dkn1=dkn1 *0.8_dl
dlnk0=dlnk0/2 !*0.3_dl
dkn2=dkn2*0.85_dl
end if
qmax_log = dkn1/dlnk0
q_switch = 2*6.3/taurst
!Want linear spacing for wavenumbers which come inside horizon
!Could use sound horizon, but for tensors that is not relevant
q_cmb = 2*l_smooth_sample/CP%chi0*AccuracyBoost !assume everything is smooth at l > l_smooth_sample
if (CP%Want_CMB .and. maximum_l > 5000 .and. CP%AccuratePolarization) q_cmb = q_cmb*1.4
!prevent EE going wild in tail
dksmooth = q_cmb/2/(AccuracyBoost)**2
if (CP%Want_CMB) dksmooth = dksmooth/6
call Ranges_Init(Evolve_q)
call Ranges_Add_delta(Evolve_q, qmin, qmax_log, dlnk0, IsLog = .true.)
call Ranges_Add_delta(Evolve_q, qmax_log, min(qmax,q_switch), dkn1)
if (qmax > q_switch) then
call Ranges_Add_delta(Evolve_q, q_switch, min(q_cmb,qmax), dkn2)
if (qmax > q_cmb) then
dksmooth = log(1 + dksmooth/q_cmb)
call Ranges_Add_delta(Evolve_q, q_cmb, qmax, dksmooth, IsLog = .true.)
end if
end if
call Ranges_GetArray(Evolve_q, .false.)
if (CP%closed) call SetClosedkValuesFromArr(Evolve_q, .false.)
end subroutine SetkValuesForSources
subroutine SetClosedkValuesFromArr(R, forInt)
Type(Regions) :: R
integer i,nu,lastnu,nmax
!nu = 3,4,5... in CP%closed case, so set nearest integers from arr array
logical, intent(in) :: forInt
integer ix
real(dl) dnu
integer, allocatable :: nu_array(:)
if (forInt .and. nint(R%points(1)*CP%r)<=3) then
!quantization is important
call Ranges_Getdpoints(R,half_ends = .false.)
R%dpoints = max(1,int(R%dpoints*CP%r+0.02))
lastnu=2
ix=1
dnu =R%dpoints(ix)