forked from cmbant/CAMB
-
Notifications
You must be signed in to change notification settings - Fork 0
/
modules.f90
3231 lines (2629 loc) · 109 KB
/
modules.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
! Modules used by cmbmain and other routines.
! Code for Anisotropies in the Microwave Background
! by Antony Lewis (http://cosmologist.info) and Anthony Challinor
! See readme.html for documentation.
!
! 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 ModelParams
use precision
use Ranges
use InitialPower
use Reionization
use Recombination
use Errors
implicit none
public
character(LEN=*), parameter :: version = 'Sept17'
integer :: FeedbackLevel = 0 !if >0 print out useful information about the model
logical :: output_file_headers = .true.
logical, parameter :: DebugMsgs=.false. !Set to true to view progress and timing
logical, parameter :: DebugEvolution = .false. !Set to true to do all the evolution for all k
real(dl) :: DebugParam = 0._dl !not used but read in, useful for parameter-dependent tests
logical :: do_bispectrum = .false.
logical, parameter :: hard_bispectrum = .false. ! e.g. warm inflation where delicate cancellations
logical, parameter :: full_bessel_integration = .false. !(go into the tails when calculating the sources)
integer, parameter :: Nu_int = 0, Nu_trunc=1, Nu_approx = 2, Nu_best = 3
!For CAMBparams%MassiveNuMethod
!Nu_int: always integrate distribution function
!Nu_trunc: switch to expansion in velocity once non-relativistic
!Nu_approx: approximate scheme - good for CMB, but not formally correct and no good for matter power
!Nu_best: automatically use mixture which is fastest and most accurate
integer, parameter :: max_Nu = 5 !Maximum number of neutrino species
integer, parameter :: max_transfer_redshifts = 150
integer, parameter :: fileio_unit = 13 !Any number not used elsewhere will do
integer, parameter :: outNone=1
integer :: max_bessels_l_index = 1000000
real(dl) :: max_bessels_etak = 1000000*2
real(dl), parameter :: OutputDenominator =twopi
!When using outNone the output is l(l+1)Cl/OutputDenominator
Type(Regions) :: TimeSteps
type TransferParams
logical :: high_precision
logical :: accurate_massive_neutrinos
integer :: num_redshifts
real(dl) :: kmax !these are acutally q values, but same as k for flat
integer :: k_per_logint ! ..
real(dl) :: redshifts(max_transfer_redshifts)
!JD 08/13 Added so both NL lensing and PK can be run at the same time
real(dl) :: PK_redshifts(max_transfer_redshifts)
real(dl) :: NLL_redshifts(max_transfer_redshifts)
integer :: PK_redshifts_index(max_transfer_redshifts)
integer :: NLL_redshifts_index(max_transfer_redshifts)
integer :: PK_num_redshifts
integer :: NLL_num_redshifts
end type TransferParams
!other variables, options, derived variables, etc.
integer, parameter :: NonLinear_none=0, NonLinear_Pk =1, NonLinear_Lens=2
integer, parameter :: NonLinear_both=3 !JD 08/13 added so both can be done
! Main parameters type
type CAMBparams
logical :: WantCls, WantTransfer
logical :: WantScalars, WantTensors, WantVectors
logical :: DoLensing
logical :: want_zstar, want_zdrag !!JH for updated BAO likelihood.
logical :: PK_WantTransfer !JD 08/13 Added so both NL lensing and PK can be run at the same time
integer :: NonLinear
logical :: Want_CMB
integer :: Max_l, Max_l_tensor
real(dl) :: Max_eta_k, Max_eta_k_tensor
! _tensor settings only used in initialization,
!Max_l and Max_eta_k are set to the tensor variables if only tensors requested
real(dl) :: omegab, omegac, omegav, omegan
!Omega baryon, CDM, Lambda and massive neutrino
real(dl) :: H0,TCMB,yhe,Num_Nu_massless
integer :: Num_Nu_massive !sum of Nu_mass_numbers below
integer :: Nu_mass_eigenstates !1 for degenerate masses
logical :: share_delta_neff !take fractional part to heat all eigenstates the same
real(dl) :: Nu_mass_degeneracies(max_nu)
real(dl) :: Nu_mass_fractions(max_nu) !The ratios of the total densities
integer :: Nu_mass_numbers(max_nu) !physical number per eigenstate
integer :: Scalar_initial_condition
!must be one of the initial_xxx values defined in GaugeInterface
integer :: OutputNormalization
!outNone, or C_OutputNormalization=1 if > 1
logical :: AccuratePolarization
!Do you care about the accuracy of the polarization Cls?
logical :: AccurateBB
!Do you care about BB accuracy (e.g. in lensing)
!Reionization settings - used if Reion%Reionization=.true.
logical :: AccurateReionization
!Do you care about pecent level accuracy on EE signal from reionization?
integer :: MassiveNuMethod
type(InitialPowerParams) :: InitPower !see power_tilt.f90 - you can change this
type(ReionizationParams) :: Reion
type(RecombinationParams):: Recomb
type(TransferParams) :: Transfer
real(dl) :: InitialConditionVector(1:10) !Allow up to 10 for future extensions
!ignored unless Scalar_initial_condition == initial_vector
logical OnlyTransfers !Don't use initial power spectrum data, instead get Delta_q_l array
!If true, sigma_8 is not calculated either
logical DerivedParameters !calculate various derived parameters (ThermoDerivedParams)
!Derived parameters, not set initially
type(ReionizationHistory) :: ReionHist
logical flat,closed,open
real(dl) omegak
real(dl) curv,r, Ksign !CP%r = 1/sqrt(|CP%curv|), CP%Ksign = 1,0 or -1
real(dl) tau0,chi0 !time today and rofChi(CP%tau0/CP%r)
end type CAMBparams
type(CAMBparams), save :: CP !Global collection of parameters
real(dl) scale !relative to CP%flat. e.g. for scaling lSamp%l sampling.
logical ::call_again = .false.
!if being called again with same parameters to get different thing
! grhom =kappa*a^2*rho_m0
! grhornomass=grhor*number of massless neutrino species
! taurst,taurend - time at start/end of recombination
! dtaurec - dtau during recombination
! adotrad - a(tau) in radiation era
real(dl) grhom,grhog,grhor,grhob,grhoc,grhov,grhornomass,grhok
real(dl) taurst,dtaurec,taurend, tau_maxvis,adotrad
!Neutrinos
real(dl) grhormass(max_nu)
! nu_masses=m_nu*c**2/(k_B*T_nu0)
real(dl) :: nu_masses(max_nu)
real(dl) akthom !sigma_T * (number density of protons now)
real(dl) fHe !n_He_tot / n_H_tot
real(dl) Nnow
integer :: ThreadNum = 0
!If zero assigned automatically, obviously only used if parallelised
!Parameters for checking/changing overall accuracy
!If HighAccuracyDefault=.false., the other parameters equal to 1 corresponds to ~0.3% scalar C_l accuracy
!If HighAccuracyDefault=.true., the other parameters equal to 1 corresponds to ~0.1% scalar C_l accuracy (at L>600)
logical :: HighAccuracyDefault = .true.
real(dl) :: lSampleBoost=1._dl
!Increase lSampleBoost to increase sampling in lSamp%l for Cl interpolation
real(dl) :: AccuracyBoost =1._dl
!Decrease step sizes, etc. by this parameter. Useful for checking accuracy.
!Can also be used to improve speed significantly if less accuracy is required.
!or improving accuracy for extreme models.
!Note this does not increase lSamp%l sampling or massive neutrino q-sampling
real(sp) :: lAccuracyBoost=1.
!Boost number of multipoles integrated in Boltzman heirarchy
integer :: limber_phiphi = 0 !for l>limber_phiphi use limber approx for lensing potential
integer :: num_redshiftwindows = 0
integer :: num_extra_redshiftwindows = 0
integer :: num_custom_sources = 0
integer, allocatable :: custom_source_ell_scales(:)
integer, parameter :: lmin = 2
!must be either 1 or 2
real(dl), parameter :: OmegaKFlat = 5e-7_dl !Value at which to use flat code
real(dl),parameter :: tol=1.0d-4 !Base tolerance for integrations
! used as parameter for spline - tells it to use 'natural' end values
real(dl), parameter :: spl_large=1.e40_dl
integer, parameter:: l0max=4000
! lmax is max possible number of l's evaluated
integer, parameter :: lmax_arr = l0max
character(LEN=1024) :: highL_unlensed_cl_template = 'HighLExtrapTemplate_lenspotentialCls.dat'
!fiducial high-accuracy high-L C_L used for making small cosmology-independent numerical corrections
!to lensing and C_L interpolation. Ideally close to models of interest, but dependence is weak.
logical :: use_spline_template = .true.
integer, parameter :: lmax_extrap_highl = 8000
real(dl), allocatable :: highL_CL_template(:,:)
integer, parameter :: derived_age=1, derived_zstar=2, derived_rstar=3, derived_thetastar=4, derived_DAstar = 5, &
derived_zdrag=6, derived_rdrag=7,derived_kD=8,derived_thetaD=9, derived_zEQ =10, derived_keq =11, &
derived_thetaEQ=12, derived_theta_rs_EQ = 13
integer, parameter :: nthermo_derived = 13
real(dl) ThermoDerivedParams(nthermo_derived)
Type TBackgroundOutputs
real(dl), pointer :: z_outputs(:) => null()
real(dl), allocatable :: H(:), DA(:), rs_by_D_v(:)
end Type TBackgroundOutputs
Type(TBackgroundOutputs), save :: BackgroundOutputs
contains
subroutine CAMBParams_Set(P, error, DoReion)
use constants
type(CAMBparams), intent(in) :: P
real(dl) GetOmegak, fractional_number, conv
integer, optional :: error !Zero if OK
logical, optional :: DoReion
logical WantReion
integer nu_i,actual_massless
real(dl) nu_massless_degeneracy, neff_i
external GetOmegak
real(dl), save :: last_tau0
!Constants in SI units
global_error_flag = 0
if ((P%WantTensors .or. P%WantVectors).and. P%WantTransfer .and. .not. P%WantScalars) then
call GlobalError( 'Cannot generate tensor C_l and transfer without scalar C_l',error_unsupported_params)
end if
if (present(error)) error = global_error_flag
if (global_error_flag/=0) return
if (present(DoReion)) then
WantReion = DoReion
else
WantReion = .true.
end if
CP=P
if (call_again) CP%DerivedParameters = .false.
CP%Max_eta_k = max(CP%Max_eta_k,CP%Max_eta_k_tensor)
if (CP%WantTransfer) then
CP%WantScalars=.true.
if (.not. CP%WantCls) then
CP%AccuratePolarization = .false.
CP%Reion%Reionization = .false.
end if
else
CP%transfer%num_redshifts=0
end if
if (CP%Num_Nu_Massive /= sum(CP%Nu_mass_numbers(1:CP%Nu_mass_eigenstates))) then
if (sum(CP%Nu_mass_numbers(1:CP%Nu_mass_eigenstates))/=0) call MpiStop('Num_Nu_Massive is not sum of Nu_mass_numbers')
end if
if (CP%Omegan == 0 .and. CP%Num_Nu_Massive /=0) then
if (CP%share_delta_neff) then
CP%Num_Nu_Massless = CP%Num_Nu_Massless + CP%Num_Nu_Massive
else
CP%Num_Nu_Massless = CP%Num_Nu_Massless + sum(CP%Nu_mass_degeneracies(1:CP%Nu_mass_eigenstates))
end if
CP%Num_Nu_Massive = 0
CP%Nu_mass_numbers = 0
end if
nu_massless_degeneracy = CP%Num_Nu_massless !N_eff for massless neutrinos
if (CP%Num_nu_massive > 0) then
if (CP%Nu_mass_eigenstates==0) call MpiStop('Have Num_nu_massive>0 but no nu_mass_eigenstates')
if (CP%Nu_mass_eigenstates==1 .and. CP%Nu_mass_numbers(1)==0) CP%Nu_mass_numbers(1) = CP%Num_Nu_Massive
if (all(CP%Nu_mass_numbers(1:CP%Nu_mass_eigenstates)==0)) CP%Nu_mass_numbers=1 !just assume one for all
if (CP%share_delta_neff) then
!default case of equal heating of all neutrinos
fractional_number = CP%Num_Nu_massless + CP%Num_Nu_massive
actual_massless = int(CP%Num_Nu_massless + 1e-6_dl)
neff_i = fractional_number/(actual_massless + CP%Num_Nu_massive)
nu_massless_degeneracy = neff_i*actual_massless
CP%Nu_mass_degeneracies(1:CP%Nu_mass_eigenstates) = CP%Nu_mass_numbers(1:CP%Nu_mass_eigenstates)*neff_i
end if
if (abs(sum(CP%Nu_mass_fractions(1:CP%Nu_mass_eigenstates))-1) > 1e-4) &
call MpiStop('Nu_mass_fractions do not add up to 1')
else
CP%Nu_mass_eigenstates = 0
end if
if ((CP%WantTransfer).and. CP%MassiveNuMethod==Nu_approx) then
CP%MassiveNuMethod = Nu_trunc
end if
CP%omegak = GetOmegak()
CP%flat = (abs(CP%omegak) <= OmegaKFlat)
CP%closed = CP%omegak < -OmegaKFlat
CP%open = .not.CP%flat.and..not.CP%closed
if (CP%flat) then
CP%curv=0
CP%Ksign=0
CP%r=1._dl !so we can use tau/CP%r, etc, where CP%r's cancel
else
CP%curv=-CP%omegak/((c/1000)/CP%h0)**2
CP%Ksign =sign(1._dl,CP%curv)
CP%r=1._dl/sqrt(abs(CP%curv))
end if
! grho gives the contribution to the expansion rate from: (g) photons,
! (r) one flavor of relativistic neutrino (2 degrees of freedom),
! (m) nonrelativistic matter (for Omega=1). grho is actually
! 8*pi*G*rho/c^2 at a=1, with units of Mpc**(-2).
! a=tau(Mpc)*adotrad, with a=1 today, assuming 3 neutrinos.
! (Used only to set the initial conformal time.)
!H0 is in km/s/Mpc
grhom = 3*CP%h0**2/c**2*1000**2 !3*h0^2/c^2 (=8*pi*G*rho_crit/c^2)
!grhom=3.3379d-11*h0*h0
grhog = kappa/c**2*4*sigma_boltz/c**3*CP%tcmb**4*Mpc**2 !8*pi*G/c^2*4*sigma_B/c^3 T^4
! grhog=1.4952d-13*tcmb**4
grhor = 7._dl/8*(4._dl/11)**(4._dl/3)*grhog !7/8*(4/11)^(4/3)*grhog (per neutrino species)
!grhor=3.3957d-14*tcmb**4
!correction for fractional number of neutrinos, e.g. 3.04 to give slightly higher T_nu hence rhor
!for massive Nu_mass_degeneracies parameters account for heating from grhor
grhornomass=grhor*nu_massless_degeneracy
grhormass=0
do nu_i = 1, CP%Nu_mass_eigenstates
grhormass(nu_i)=grhor*CP%Nu_mass_degeneracies(nu_i)
end do
grhoc=grhom*CP%omegac
grhob=grhom*CP%omegab
grhov=grhom*CP%omegav
grhok=grhom*CP%omegak
! adotrad gives the relation a(tau) in the radiation era:
adotrad = sqrt((grhog+grhornomass+sum(grhormass(1:CP%Nu_mass_eigenstates)))/3)
Nnow = CP%omegab*(1-CP%yhe)*grhom*c**2/kappa/m_H/Mpc**2
akthom = sigma_thomson*Nnow*Mpc
!sigma_T * (number density of protons now)
fHe = CP%YHe/(mass_ratio_He_H*(1.d0-CP%YHe)) !n_He_tot / n_H_tot
if (.not.call_again) then
call init_massive_nu(CP%omegan /=0)
call init_background
if (global_error_flag==0) then
CP%tau0=TimeOfz(0._dl)
! print *, 'chi = ', (CP%tau0 - TimeOfz(0.15_dl)) * CP%h0/100
last_tau0=CP%tau0
if (WantReion) call Reionization_Init(CP%Reion,CP%ReionHist, CP%YHe, akthom, CP%tau0, FeedbackLevel)
end if
else
CP%tau0=last_tau0
end if
!JD 08/13 Changes for nonlinear lensing of CMB + MPK compatibility
!if ( CP%NonLinear==NonLinear_Lens) then
if (CP%NonLinear==NonLinear_Lens .or. CP%NonLinear==NonLinear_both ) then
CP%Transfer%kmax = max(CP%Transfer%kmax, CP%Max_eta_k/CP%tau0)
if (FeedbackLevel > 0 .and. CP%Transfer%kmax== CP%Max_eta_k/CP%tau0) &
write (*,*) 'max_eta_k changed to ', CP%Max_eta_k
end if
if (CP%closed .and. CP%tau0/CP%r >3.14) then
call GlobalError('chi >= pi in closed model not supported',error_unsupported_params)
end if
if (global_error_flag/=0) then
if (present(error)) error = global_error_flag
return
end if
if (present(error)) then
error = 0
else if (FeedbackLevel > 0 .and. .not. call_again) then
write(*,'("Om_b h^2 = ",f9.6)') CP%omegab*(CP%H0/100)**2
write(*,'("Om_c h^2 = ",f9.6)') CP%omegac*(CP%H0/100)**2
write(*,'("Om_nu h^2 = ",f9.6)') CP%omegan*(CP%H0/100)**2
write(*,'("Om_Lambda = ",f9.6)') CP%omegav
write(*,'("Om_K = ",f9.6)') CP%omegak
write(*,'("Om_m (1-Om_K-Om_L) = ",f9.6)') 1-CP%omegak-CP%omegav
write(*,'("100 theta (CosmoMC) = ",f9.6)') 100*CosmomcTheta()
if (CP%Num_Nu_Massive > 0) then
write(*,'("N_eff (total) = ",f9.6)') nu_massless_degeneracy + &
sum(CP%Nu_mass_degeneracies(1:CP%Nu_mass_eigenstates))
do nu_i=1, CP%Nu_mass_eigenstates
conv = k_B*(8*grhor/grhog/7)**0.25*CP%tcmb/eV * &
(CP%nu_mass_degeneracies(nu_i)/CP%nu_mass_numbers(nu_i))**0.25 !approx 1.68e-4
write(*,'(I2, " nu, g=",f7.4," m_nu*c^2/k_B/T_nu0= ",f9.2," (m_nu= ",f6.3," eV)")') &
CP%nu_mass_numbers(nu_i), CP%nu_mass_degeneracies(nu_i), nu_masses(nu_i),conv*nu_masses(nu_i)
end do
end if
end if
CP%chi0=rofChi(CP%tau0/CP%r)
scale= CP%chi0*CP%r/CP%tau0 !e.g. change l sampling depending on approx peak spacing
end subroutine CAMBParams_Set
function GetTestTime()
real(sp) GetTestTime
real(sp) atime
! GetTestTime = etime(tarray)
!Can replace this if etime gives problems
!Or just comment out - only used if DebugMsgs = .true.
call cpu_time(atime)
GetTestTime = atime
end function GetTestTime
function rofChi(Chi) !sinh(chi) for open, sin(chi) for closed.
real(dl) Chi,rofChi
if (CP%closed) then
rofChi=sin(chi)
else if (CP%open) then
rofChi=sinh(chi)
else
rofChi=chi
endif
end function rofChi
function cosfunc (Chi)
real(dl) Chi,cosfunc
if (CP%closed) then
cosfunc= cos(chi)
else if (CP%open) then
cosfunc=cosh(chi)
else
cosfunc = 1._dl
endif
end function cosfunc
function tanfunc(Chi)
real(dl) Chi,tanfunc
if (CP%closed) then
tanfunc=tan(Chi)
else if (CP%open) then
tanfunc=tanh(Chi)
else
tanfunc=Chi
end if
end function tanfunc
function invsinfunc(x)
real(dl) invsinfunc,x
if (CP%closed) then
invsinfunc=asin(x)
else if (CP%open) then
invsinfunc=log((x+sqrt(1._dl+x**2)))
else
invsinfunc = x
endif
end function invsinfunc
function f_K(x)
real(dl) :: f_K
real(dl), intent(in) :: x
f_K = CP%r*rofChi(x/CP%r)
end function f_K
function DeltaTime(a1,a2, in_tol)
implicit none
real(dl) DeltaTime, atol
real(dl), intent(IN) :: a1,a2
real(dl), optional, intent(in) :: in_tol
real(dl) dtauda, rombint !diff of tau w.r.t a and integration
external dtauda, rombint
if (present(in_tol)) then
atol = in_tol
else
atol = tol/1000/exp(AccuracyBoost-1)
end if
DeltaTime=rombint(dtauda,a1,a2,atol)
end function DeltaTime
function TimeOfz(z)
implicit none
real(dl) TimeOfz
real(dl), intent(IN) :: z
TimeOfz=DeltaTime(0._dl,1._dl/(z+1._dl))
end function TimeOfz
subroutine TimeOfzArr(nz, redshifts, outputs)
integer, intent(in) :: nz
real(dl), intent(in) :: redshifts(nz)
real(dl), intent(out) :: outputs(nz)
integer i
!Dumb slow version
!$OMP PARALLEL DO DEFAUlT(SHARED)
do i=1, nz
outputs(i) = timeOfZ(redshifts(i))
end do
!$OMP END PARALLEL DO
end subroutine TimeOfzArr
function DeltaPhysicalTimeGyr(a1,a2, in_tol)
use constants
real(dl), intent(in) :: a1, a2
real(dl), optional, intent(in) :: in_tol
real(dl) rombint,DeltaPhysicalTimeGyr, atol
external rombint
if (present(in_tol)) then
atol = in_tol
else
atol = 1d-4/exp(AccuracyBoost-1)
end if
DeltaPhysicalTimeGyr = rombint(dtda,a1,a2,atol)*Mpc/c/Gyr
end function DeltaPhysicalTimeGyr
function AngularDiameterDistance(z)
!This is the physical (non-comoving) angular diameter distance in Mpc
real(dl) AngularDiameterDistance
real(dl), intent(in) :: z
AngularDiameterDistance = CP%r/(1+z)*rofchi(ComovingRadialDistance(z) /CP%r)
end function AngularDiameterDistance
subroutine AngularDiameterDistanceArr(arr, z, n)
!This is the physical (non-comoving) angular diameter distance in Mpc for array of z
!z array must be monotonically increasing
integer,intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: z(n)
integer i
call ComovingRadialDistanceArr(arr, z, n, 1e-4_dl)
if (CP%flat) then
arr = arr/(1+z)
else
do i=1, n
arr(i) = CP%r/(1+z(i))*rofchi(arr(i)/CP%r)
end do
end if
end subroutine AngularDiameterDistanceArr
function AngularDiameterDistance2(z1, z2) ! z1 < z2
!From http://www.slac.stanford.edu/~amantz/work/fgas14/#cosmomc
real(dl) AngularDiameterDistance2
real(dl), intent(in) :: z1, z2
AngularDiameterDistance2 = CP%r/(1+z2)*rofchi(ComovingRadialDistance(z2)/CP%r - ComovingRadialDistance(z1)/CP%r)
end function AngularDiameterDistance2
function LuminosityDistance(z)
real(dl) LuminosityDistance
real(dl), intent(in) :: z
LuminosityDistance = AngularDiameterDistance(z)*(1+z)**2
end function LuminosityDistance
function ComovingRadialDistance(z)
real(dl) ComovingRadialDistance
real(dl), intent(in) :: z
ComovingRadialDistance = DeltaTime(1/(1+z),1._dl)
end function ComovingRadialDistance
subroutine ComovingRadialDistanceArr(arr, z, n, tol)
!z array must be monotonically increasing
integer, intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: z(n)
real(dl), intent(in) :: tol
integer i
!$OMP PARALLEL DO DEFAULT(SHARED),SCHEDULE(STATIC)
do i = 1, n
if (i==1) then
if (z(i) < 1e-6_dl) then
arr(i) = 0
else
arr(i) = DeltaTime(1/(1+z(i)),1._dl, tol)
end if
else
if (z(i) < z(i-1)) error stop 'ComovingRadialDistanceArr redshifts out of order'
arr(i) = DeltaTime(1/(1+z(i)),1/(1+z(i-1)),tol)
end if
end do
!$OMP END PARALLEL DO
do i = 2, n
arr(i) = arr(i) + arr(i-1)
end do
end subroutine ComovingRadialDistanceArr
function Hofz(z)
!non-comoving Hubble in MPC units, divide by MPC_in_sec to get in SI units
!multiply by c/1e3 to get in km/s/Mpc units
real(dl) Hofz, dtauda,a
real(dl), intent(in) :: z
external dtauda
a = 1/(1+z)
Hofz = 1/(a**2*dtauda(a))
end function Hofz
subroutine HofzArr(arr, z, n)
!non-comoving Hubble in MPC units, divide by MPC_in_sec to get in SI units
!multiply by c/1e3 to get in km/s/Mpc units
integer,intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: z(n)
integer i
do i=1, n
arr(i) = Hofz(z(i))
end do
end subroutine HofzArr
real(dl) function BAO_D_v_from_DA_H(z, DA, Hz)
real(dl), intent(in) :: z, DA, Hz
real(dl) ADD
ADD = DA*(1.d0+z)
BAO_D_v_from_DA_H = ((ADD)**2.d0*z/Hz)**(1.d0/3.d0)
end function BAO_D_v_from_DA_H
real(dl) function BAO_D_v(z)
real(dl), intent(IN) :: z
BAO_D_v = BAO_D_v_from_DA_H(z,AngularDiameterDistance(z), Hofz(z))
end function BAO_D_v
function dsound_da_exact(a)
implicit none
real(dl) dsound_da_exact,dtauda,a,R,cs
external dtauda
R = 3*grhob*a / (4*grhog)
cs=1.0d0/sqrt(3*(1+R))
dsound_da_exact=dtauda(a)*cs
end function dsound_da_exact
function dsound_da(a)
!approximate form used e.g. by CosmoMC for theta
implicit none
real(dl) dsound_da,dtauda,a,R,cs
external dtauda
R=3.0d4*a*CP%omegab*(CP%h0/100.0d0)**2
! R = 3*grhob*a / (4*grhog) //above is mostly within 0.2% and used for previous consistency
cs=1.0d0/sqrt(3*(1+R))
dsound_da=dtauda(a)*cs
end function dsound_da
function dtda(a)
real(dl) dtda,dtauda,a
external dtauda
dtda= dtauda(a)*a
end function
function CosmomcTheta()
real(dl) zstar, astar, atol, rs, DA
real(dl) CosmomcTheta
real(dl) ombh2, omdmh2
real(dl) rombint
external rombint
ombh2 = CP%omegab*(CP%h0/100.0d0)**2
omdmh2 = (CP%omegac+CP%omegan)*(CP%h0/100.0d0)**2
!!From Hu & Sugiyama
zstar = 1048*(1+0.00124*ombh2**(-0.738))*(1+ &
(0.0783*ombh2**(-0.238)/(1+39.5*ombh2**0.763)) * &
(omdmh2+ombh2)**(0.560/(1+21.1*ombh2**1.81)))
astar = 1/(1+zstar)
atol = 1e-6
rs = rombint(dsound_da,1d-8,astar,atol)
DA = AngularDiameterDistance(zstar)/astar
CosmomcTheta = rs/DA
! print *,'z* = ',zstar, 'r_s = ',rs, 'DA = ',DA, rs/DA
end function CosmomcTheta
end module ModelParams
!ccccccccccccccccccccccccccccccccccccccccccccccccccc
module lvalues
use precision
use ModelParams
implicit none
public
Type lSamples
integer l0
integer l(lmax_arr)
end Type lSamples
Type(lSamples) :: lSamp
!Sources
logical :: Log_lvalues = .false.
contains
function lvalues_indexOf(lSet,l)
type(lSamples) :: lSet
integer, intent(in) :: l
integer lvalues_indexOf, i
do i=2,lSet%l0
if (l < lSet%l(i)) then
lvalues_indexOf = i-1
return
end if
end do
lvalues_indexOf = lSet%l0
end function lvalues_indexOf
subroutine initlval(lSet,max_l)
! This subroutines initializes lSet%l arrays. Other values will be interpolated.
implicit none
type(lSamples) :: lSet
integer, intent(IN) :: max_l
integer lind, lvar, step,top,bot,ls(lmax_arr)
real(dl) AScale
Ascale=scale/lSampleBoost
if (lSampleBoost >=50) then
!just do all of them
lind=0
do lvar=lmin, max_l
lind=lind+1
ls(lind)=lvar
end do
lSet%l0=lind
lSet%l(1:lind) = ls(1:lind)
return
end if
lind=0
do lvar=lmin, 10
lind=lind+1
ls(lind)=lvar
end do
if (CP%AccurateReionization) then
if (lSampleBoost > 1) then
do lvar=11, 37,1
lind=lind+1
ls(lind)=lvar
end do
else
do lvar=11, 37,2
lind=lind+1
ls(lind)=lvar
end do
end if
step = max(nint(5*Ascale),2)
bot=40
top=bot + step*10
else
if (lSampleBoost >1) then
do lvar=11, 15
lind=lind+1
ls(lind)=lvar
end do
else
lind=lind+1
ls(lind)=12
lind=lind+1
ls(lind)=15
end if
step = max(nint(10*Ascale),3)
bot=15+max(step/2,2)
top=bot + step*7
end if
do lvar=bot, top, step
lind=lind+1
ls(lind)=lvar
end do
!Sources
if (Log_lvalues) then
!Useful for generating smooth things like 21cm to high l
step=max(nint(20*Ascale),4)
do
lvar = lvar + step
if (lvar > max_l) exit
lind=lind+1
ls(lind)=lvar
step = nint(step*1.2) !log spacing
end do
else
step=max(nint(20*Ascale),4)
bot=ls(lind)+step
top=bot+step*2
do lvar = bot,top,step
lind=lind+1
ls(lind)=lvar
end do
if (ls(lind)>=max_l) then
do lvar=lind,1,-1
if (ls(lvar)<=max_l) exit
end do
lind=lvar
if (ls(lind)<max_l) then
lind=lind+1
ls(lind)=max_l
end if
else
step=max(nint(25*Ascale),4)
!Get EE right around l=200 by putting extra point at 175
bot=ls(lind)+step
top=bot+step
do lvar = bot,top,step
lind=lind+1
ls(lind)=lvar
end do
if (ls(lind)>=max_l) then
do lvar=lind,1,-1
if (ls(lvar)<=max_l) exit
end do
lind=lvar
if (ls(lind)<max_l) then
lind=lind+1
ls(lind)=max_l
end if
else
if (HighAccuracyDefault .and. .not. use_spline_template) then
step=max(nint(42*Ascale),7)
else
step=max(nint(50*Ascale),7)
end if
bot=ls(lind)+step
top=min(5000,max_l)
do lvar = bot,top,step
lind=lind+1
ls(lind)=lvar
end do
if (max_l > 5000) then
!Should be pretty smooth or tiny out here
step=max(nint(400*Ascale),50)
lvar = ls(lind)
do
lvar = lvar + step
if (lvar > max_l) exit
lind=lind+1
ls(lind)=lvar
step = nint(step*1.5) !log spacing
end do
end if
!Sources
end if !log_lvalues
if (ls(lind) /=max_l) then
lind=lind+1
ls(lind)=max_l
end if
if (.not. CP%flat) ls(lind-1)=int(max_l+ls(lind-2))/2
!Not in CP%flat case so interpolation table is the same when using lower l_max
end if
end if
lSet%l0=lind
lSet%l(1:lind) = ls(1:lind)
end subroutine initlval
subroutine InterpolateClArr(lSet,iCl, all_Cl, max_ind)
type (lSamples), intent(in) :: lSet
real(dl), intent(in) :: iCl(*)
real(dl), intent(out):: all_Cl(lmin:*)
integer, intent(in) :: max_ind
integer il,llo,lhi, xi
real(dl) ddCl(lSet%l0)
real(dl) xl(lSet%l0)
real(dl) a0,b0,ho
real(dl), parameter :: cllo=1.e30_dl,clhi=1.e30_dl
if (max_ind > lSet%l0) call MpiStop('Wrong max_ind in InterpolateClArr')
xl = real(lSet%l(1:lSet%l0),dl)
call spline(xl,iCL(1),max_ind,cllo,clhi,ddCl(1))
llo=1
do il=lmin,lSet%l(max_ind)
xi=il
if ((xi > lSet%l(llo+1)).and.(llo < max_ind)) then
llo=llo+1
end if
lhi=llo+1
ho=lSet%l(lhi)-lSet%l(llo)
a0=(lSet%l(lhi)-xi)/ho
b0=(xi-lSet%l(llo))/ho
all_Cl(il) = a0*iCl(llo)+ b0*iCl(lhi)+((a0**3-a0)* ddCl(llo) &
+(b0**3-b0)*ddCl(lhi))*ho**2/6
end do
end subroutine InterpolateClArr
subroutine InterpolateClArrTemplated(lSet,iCl, all_Cl, max_ind, template_index)
type (lSamples), intent(in) :: lSet
real(dl), intent(in) :: iCl(*)
real(dl), intent(out):: all_Cl(lmin:*)
integer, intent(in) :: max_ind
integer, intent(in), optional :: template_index
integer maxdelta, il
real(dl) DeltaCL(lSet%l0)
real(dl), allocatable :: tmpall(:)
if (max_ind > lSet%l0) call MpiStop('Wrong max_ind in InterpolateClArrTemplated')
if (use_spline_template .and. present(template_index)) then
if (template_index<=3) then
!interpolate only the difference between the C_l and an accurately interpolated template. Temp only for the mo.