-
Notifications
You must be signed in to change notification settings - Fork 159
/
Copy pathresults.f90
4107 lines (3479 loc) · 156 KB
/
results.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 results
use constants, only : const_pi, const_twopi
use MiscUtils
use RangeUtils
use StringUtils
use MathUtils
use config
use model
use splines
implicit none
public
Type TBackgroundOutputs
real(dl), allocatable :: H(:), DA(:), rs_by_D_v(:)
end Type TBackgroundOutputs
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
Type lSamples
integer :: nl = 0
integer :: lmin = 2
integer, allocatable :: l(:)
logical :: use_spline_template = .true.
contains
procedure :: Init => lSamples_init
procedure :: IndexOf => lSamples_indexOf
procedure :: InterpolateClArr
procedure :: InterpolateClArrTemplated
end Type lSamples
logical, parameter :: dowinlens = .false. !not used, test getting CMB lensing using visibility
integer, parameter :: thermal_history_def_timesteps = 20000
Type TThermoData
logical :: HasThermoData = .false. !Has it been computed yet for current parameters?
!Background thermal history, interpolated from precomputed tables
integer :: nthermo !Number of table steps
!baryon temperature, sound speed, ionization fractions, and opacity
real(dl), dimension(:), allocatable :: tb, cs2, xe, dotmu
! e^(-tau) and derivatives
real(dl), dimension(:), allocatable :: emmu, dcs2,demmu, ddotmu, dddotmu, ddddotmu
real(dl), dimension(:), allocatable :: ScaleFactor, dScaleFactor, adot, dadot
real(dl), dimension(:), allocatable :: winlens, dwinlens
real(dl) tauminn,dlntau
real(dl) :: tight_tau, actual_opt_depth
!Times when 1/(opacity*tau) = 0.01, for use switching tight coupling approximation
real(dl) :: matter_verydom_tau
real(dl) :: recombination_saha_tau
!sound horizon and recombination redshifts
real(dl) :: r_drag0, z_star, z_drag !!JH for updated BAO likelihood.
real(dl), dimension(:), allocatable :: step_redshift, rhos_fac, drhos_fac
real(dl) :: tau_start_redshiftwindows,tau_end_redshiftwindows
logical :: has_lensing_windows = .false.
real(dl) recombination_Tgas_tau
Type(TCubicSpline) :: ScaleFactorAtTime
!Mapping between redshift and time
real(dl), private,dimension(:), allocatable :: redshift_time, dredshift_time
real(dl), private, dimension(:), allocatable :: arhos_fac, darhos_fac, ddarhos_fac
contains
procedure :: Init => Thermo_Init
procedure :: OpacityToTime => Thermo_OpacityToTime
procedure :: values => Thermo_values
procedure :: expansion_values => Thermo_expansion_values
procedure :: IonizationFunctionsAtTime
procedure, private :: DoWindowSpline
procedure, private :: SetTimeSteps
procedure, private :: SetTimeStepWindows
end type TThermoData
!Sources
Type CalWins
real(dl), allocatable :: awin_lens(:), dawin_lens(:)
end Type CalWins
Type LimberRec
integer n1,n2 !corresponding time step array indices
real(dl), dimension(:), allocatable :: k
real(dl), dimension(:), allocatable :: Source
end Type LimberRec
Type ClTransferData
!Cl transfer function variables
!values of q for integration over q to get C_ls
Type (lSamples) :: ls ! l that are computed
integer :: NumSources
!Changes -scalars: 2 for just CMB, 3 for lensing
!- tensors: T and E and phi (for lensing), and T, E, B respectively
type (TRanges) :: q
real(dl), dimension(:,:,:), allocatable :: Delta_p_l_k
!The L index of the lowest L to use for Limber
integer, dimension(:), allocatable :: Limber_l_min
!For each l, the set of k in each limber window
!indices LimberWindow(SourceNum,l)
Type(LimberRec), dimension(:,:), allocatable :: Limber_windows
!The maximum L needed for non-Limber
integer max_index_nonlimber
end Type ClTransferData
type TCLdata
Type(ClTransferData) :: CTransScal, CTransTens, CTransVec
real(dl), dimension (:,:), allocatable :: Cl_scalar, Cl_tensor, Cl_vector
!Indices are Cl_xxx( l , Cl_type)
!where Cl_type is one of the above constants
real(dl), dimension (:,:,:), allocatable :: Cl_Scalar_Array
!Indices are Cl_xxx( l , field1,field2)
!where ordering of fields is T, E, \psi (CMB lensing potential), window_1, window_2...
!The following are set only if doing lensing
integer lmax_lensed !Only accurate to rather less than this
real(dl) , dimension (:,:), allocatable :: Cl_lensed
!Cl_lensed(l, Cl_type) are the interpolated Cls
contains
procedure :: InitCls => TCLdata_InitCls
procedure :: output_cl_files => TCLdata_output_cl_files
procedure :: output_lens_pot_files => TCLdata_output_lens_pot_files
procedure :: NormalizeClsAtL => TCLdata_NormalizeClsAtL
procedure :: output_veccl_files => TCLdata_output_veccl_files
end type TCLdata
Type TTimeSources
! values of q to evolve the propagation equations to compute the sources
type(TRanges) :: Evolve_q
real(dl), dimension(:,:,:), allocatable :: LinearSrc !Sources and second derivs
!LinearSrc indices ( k_index, source_index, time_step_index )
integer SourceNum, NonCustomSourceNum
!SourceNum is total number sources (2 or 3 for scalars, 3 for tensors).
end type TTimeSources
type, extends(TCAMBdata) :: CAMBdata
type(CAMBparams) :: CP
real(dl) ThermoDerivedParams(nthermo_derived)
logical flat,closed
! grhocrit =kappa*a^2*rho_crit(0)
! 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) grhocrit,grhog,grhor,grhob,grhoc,grhov,grhornomass,grhok
real(dl) taurst,dtaurec,taurend,tau_maxvis,adotrad
real(dl) Omega_de
real(dl) curv, curvature_radius, Ksign !curvature_radius = 1/sqrt(|curv|), Ksign = 1,0 or -1
real(dl) tau0,chi0 !time today and rofChi(tau0/curvature_radius)
real(dl) scale !relative to flat. e.g. for scaling lSamp%l sampling.
real(dl) akthom !sigma_T * (number density of protons now)
real(dl) fHe !n_He_tot / n_H_tot
real(dl) Nnow
real(dl) z_eq !defined assuming all neutrinos massless
!Neutrinos
real(dl) grhormass(max_nu)
! nu_masses=m_nu*c**2/(k_B*T_nu0)
real(dl) nu_masses(max_nu)
integer :: num_transfer_redshifts = 1
real(dl), allocatable :: transfer_redshifts(:)
integer :: PK_redshifts_index(max_transfer_redshifts)
logical :: OnlyTransfer = .false. !C_L/PK not computed; initial power spectrum data, instead get Delta_q_l array
!If true, sigma_8 is not calculated either]]
logical :: HasScalarTimeSources = .false. !No power spectra, only time transfer functions
logical :: get_growth_sigma8 = .true.
!gets sigma_vdelta, like sigma8 but using velocity-density cross power,
!in late LCDM f*sigma8 = sigma_vdelta^2/sigma8
logical :: needs_good_pk_sampling = .false.
logical ::call_again = .false.
!if being called again with same parameters to get different thing
real(dl) :: reion_tau_start, reion_tau_complete
integer :: reion_n_steps
Type(TNuPerturbations) :: NuPerturbations
Type(TBackgroundOutputs) :: BackgroundOutputs
!Time steps for sampling sources
Type(TRanges) :: TimeSteps
!Background interpolation tables for thermal history etc.
Type(TThermoData) :: ThermoData
real(dl), allocatable :: transfer_times(:)
!Matter transfer data
Type (MatterTransferData):: MT
!Matter power spectrum for default variable (used for non-linear corrections)
Type(MatterPowerData), allocatable :: CAMB_PK
Type(TClData) :: CLdata
integer :: num_redshiftwindows = 0
integer :: num_extra_redshiftwindows = 0
Type(TRedWin), allocatable :: Redshift_W(:)
real(dl), dimension(:), allocatable :: optical_depths_for21cm
Type(TTimeSources), allocatable :: ScalarTimeSources
integer :: Scalar_C_last = C_PhiE
contains
procedure :: DeltaTime => CAMBdata_DeltaTime
procedure :: DeltaTimeArr => CAMBdata_DeltaTimeArr
procedure :: TimeOfz => CAMBdata_TimeOfz
procedure :: TimeOfzArr => CAMBdata_TimeOfzArr
procedure :: DeltaPhysicalTimeGyr => CAMBdata_DeltaPhysicalTimeGyr
procedure :: DeltaPhysicalTimeGyrArr => CAMBdata_DeltaPhysicalTimeGyrArr
procedure :: AngularDiameterDistance => CAMBdata_AngularDiameterDistance
procedure :: AngularDiameterDistanceArr => CAMBdata_AngularDiameterDistanceArr
procedure :: AngularDiameterDistance2 => CAMBdata_AngularDiameterDistance2
procedure :: AngularDiameterDistance2Arr => CAMBdata_AngularDiameterDistance2Arr
procedure :: LuminosityDistance => CAMBdata_LuminosityDistance
procedure :: ComovingRadialDistance => CAMBdata_ComovingRadialDistance
procedure :: ComovingRadialDistanceArr => CAMBdata_ComovingRadialDistanceArr
procedure :: GetBackgroundDensities => CAMBdata_GetBackgroundDensities
procedure :: Hofz => CAMBdata_Hofz
procedure :: HofzArr => CAMBdata_HofzArr
procedure :: sound_horizon => CAMBdata_sound_horizon
procedure :: sound_horizon_zArr => CAMBdata_sound_horizon_zArr
procedure :: RedshiftAtTimeArr => CAMBdata_RedshiftAtTimeArr
procedure :: BAO_D_v => CAMBdata_BAO_D_v
procedure :: CosmomcTheta => CAMBdata_CosmomcTheta
procedure :: get_lmax_lensed => CAMBdata_get_lmax_lensed
procedure :: get_zstar => CAMBdata_get_zstar
procedure :: DarkEnergyStressEnergy => CAMBdata_DarkEnergyStressEnergy
procedure :: SetParams => CAMBdata_SetParams
procedure :: Free => CAMBdata_Free
procedure :: grho_no_de
procedure :: GetReionizationOptDepth
procedure :: rofChi
procedure :: cosfunc
procedure :: tanfunc
procedure :: invsinfunc
procedure :: GetComputedPKRedshifts
procedure :: binary_search
procedure, nopass :: PythonClass => CAMBdata_PythonClass
procedure, nopass :: SelfPointer => CAMBdata_SelfPointer
end type CAMBdata
interface
FUNCTION state_function(obj, a)
use precision
import
class(CAMBdata) :: obj
real(dl), intent(in) :: a
real(dl) :: state_function
END FUNCTION state_function
end interface
procedure(obj_function), private :: dtauda
contains
function CAMBdata_PythonClass()
character(LEN=:), allocatable :: CAMBdata_PythonClass
CAMBdata_PythonClass = 'CAMBdata'
end function CAMBdata_PythonClass
subroutine CAMBdata_SelfPointer(cptr,P)
use iso_c_binding
Type(c_ptr) :: cptr
Type (CAMBdata), pointer :: PType
class (TPythonInterfacedClass), pointer :: P
call c_f_pointer(cptr, PType)
P => PType
end subroutine CAMBdata_SelfPointer
subroutine CAMBdata_SetParams(this, P, error, DoReion, call_again, background_only)
!Initialize background variables; does not yet calculate thermal history
use constants
class(CAMBdata), target :: this
type(CAMBparams), intent(in) :: P
real(dl) fractional_number, conv
integer, optional :: error !Zero if OK
logical, optional :: DoReion
logical, optional :: call_again, background_only
logical WantReion, calling_again
integer nu_i,actual_massless
real(dl) nu_massless_degeneracy, neff_i, eta_k, h2
real(dl) zpeak, sigma_z, zpeakstart, zpeakend
Type(TRedWin), pointer :: Win
logical back_only
!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 (.not. allocated(P%DarkEnergy)) then
call GlobalError('DarkEnergy not set', error_darkenergy)
end if
if (present(error)) error = global_error_flag
if (global_error_flag/=0) return
WantReion = DefaultTrue(DoReion)
calling_again= DefaultFalse(call_again)
back_only = DefaultFalse(background_only)
if (calling_again) then
this%CP%WantDerivedParameters = .false.
this%CP%Accuracy = P%Accuracy
this%CP%Transfer%high_precision = P%Transfer%high_precision
this%CP%Reion%Reionization = P%Reion%Reionization
this%CP%WantTransfer =P%WantTransfer
this%CP%WantScalars =P%WantScalars
this%CP%WantTensors =P%WantTensors
this%CP%WantVectors =P%WantVectors
this%CP%WantCls = P%WantCls
else
this%CP=P
this%CP%Max_eta_k = max(this%CP%Max_eta_k,this%CP%Max_eta_k_tensor)
end if
if (P%WantTransfer .and. .not. back_only) then
this%CP%WantScalars=.true.
if (.not. P%WantCls) then
this%CP%Accuracy%AccuratePolarization = .false.
!Sources
this%CP%Reion%Reionization = this%CP%transfer_21cm_cl
end if
call this%GetComputedPKRedshifts(this%CP)
end if
if (this%CP%WantTransfer.and. this%CP%MassiveNuMethod==Nu_approx) then
this%CP%MassiveNuMethod = Nu_trunc
end if
if (.not. calling_again) then
this%ThermoData%HasThermoData = .false.
if (this%CP%Num_Nu_Massive /= sum(this%CP%Nu_mass_numbers(1:this%CP%Nu_mass_eigenstates))) then
if (sum(this%CP%Nu_mass_numbers(1:this%CP%Nu_mass_eigenstates))/=0) &
call GlobalError('Num_Nu_Massive is not sum of Nu_mass_numbers', error_unsupported_params)
end if
10 if (this%CP%Omnuh2 < 1.e-7_dl) this%CP%Omnuh2 = 0
if (this%CP%Omnuh2==0 .and. this%CP%Num_Nu_Massive /=0) then
if (this%CP%share_delta_neff) then
this%CP%Num_Nu_Massless = this%CP%Num_Nu_Massless + this%CP%Num_Nu_Massive
else
this%CP%Num_Nu_Massless = this%CP%Num_Nu_Massless + sum(this%CP%Nu_mass_degeneracies(1:this%CP%Nu_mass_eigenstates))
end if
this%CP%Num_Nu_Massive = 0
this%CP%Nu_mass_numbers = 0
end if
nu_massless_degeneracy = this%CP%Num_Nu_massless !N_eff for massless neutrinos
if (this%CP%Num_nu_massive > 0) then
if (this%CP%Nu_mass_eigenstates==0) &
call GlobalError('Have Num_nu_massive>0 but no nu_mass_eigenstates', error_unsupported_params)
if (this%CP%Nu_mass_eigenstates==1 .and. this%CP%Nu_mass_numbers(1)==0) &
this%CP%Nu_mass_numbers(1) = this%CP%Num_Nu_Massive
if (all(this%CP%Nu_mass_numbers(1:this%CP%Nu_mass_eigenstates)==0)) this%CP%Nu_mass_numbers=1 !just assume one for all
if (this%CP%share_delta_neff) then
!default case of equal heating of all neutrinos
fractional_number = this%CP%Num_Nu_massless + this%CP%Num_Nu_massive
actual_massless = int(this%CP%Num_Nu_massless + 1e-6_dl)
neff_i = fractional_number/(actual_massless + this%CP%Num_Nu_massive)
nu_massless_degeneracy = neff_i*actual_massless
this%CP%Nu_mass_degeneracies(1:this%CP%Nu_mass_eigenstates) = &
this%CP%Nu_mass_numbers(1:this%CP%Nu_mass_eigenstates)*neff_i
end if
if (abs(sum(this%CP%Nu_mass_fractions(1:this%CP%Nu_mass_eigenstates))-1) > 1e-4) &
call GlobalError('Nu_mass_fractions do not add up to 1', error_unsupported_params)
else
this%CP%Nu_mass_eigenstates = 0
end if
if (global_error_flag/=0) then
if (present(error)) error = global_error_flag
return
end if
call This%ThermoData%ScaleFactorAtTime%Clear()
this%flat = (abs(this%CP%omk) <= OmegaKFlat)
this%closed = this%CP%omk < -OmegaKFlat
if (this%flat) then
this%curv=0
this%Ksign=0
this%curvature_radius=1._dl !so we can use tau/curvature_radius, etc, where r's cancel
else
this%curv=-this%CP%omk/((c/1000)/this%CP%h0)**2
this%Ksign =sign(1._dl,this%curv)
this%curvature_radius=1._dl/sqrt(abs(this%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
this%grhocrit = 3*this%CP%h0**2/c**2*1000**2 !3*h0^2/c^2 (=8*pi*G*rho_crit/c^2)
this%grhog = kappa/c**2*4*sigma_boltz/c**3*this%CP%tcmb**4*Mpc**2 !8*pi*G/c^2*4*sigma_B/c^3 T^4
! grhog=1.4952d-13*tcmb**4
this%grhor = 7._dl/8*(4._dl/11)**(4._dl/3)*this%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
this%grhornomass=this%grhor*nu_massless_degeneracy
this%grhormass=0
do nu_i = 1, this%CP%Nu_mass_eigenstates
this%grhormass(nu_i)=this%grhor*this%CP%Nu_mass_degeneracies(nu_i)
end do
h2 = (this%CP%H0/100)**2
this%grhoc=this%grhocrit*this%CP%omch2/h2
this%grhob=this%grhocrit*this%CP%ombh2/h2
this%grhok=this%grhocrit*this%CP%omk
this%Omega_de = 1 -(this%CP%omch2 + this%CP%ombh2 + this%CP%omnuh2)/h2 - this%CP%omk &
- (this%grhornomass + this%grhog)/this%grhocrit
this%grhov=this%grhocrit*this%Omega_de
! adotrad gives da/dtau in the asymptotic radiation-dominated era:
this%adotrad = sqrt((this%grhog+this%grhornomass+sum(this%grhormass(1:this%CP%Nu_mass_eigenstates)))/3)
this%Nnow = this%CP%ombh2/h2*(1-this%CP%yhe)*this%grhocrit*c**2/kappa/m_H/Mpc**2
this%akthom = sigma_thomson*this%Nnow*Mpc
!sigma_T * (number density of protons now)
this%fHe = this%CP%YHe/(mass_ratio_He_H*(1.d0-this%CP%YHe)) !n_He_tot / n_H_tot
this%z_eq = (this%grhob+this%grhoc)/(this%grhog+this%grhornomass+sum(this%grhormass(1:this%CP%Nu_mass_eigenstates))) -1
if (this%CP%omnuh2/=0) then
!Initialize things for massive neutrinos
call ThermalNuBackground%Init()
call this%NuPerturbations%Init(P%Accuracy%AccuracyBoost*P%Accuracy%neutrino_q_boost)
! nu_masses=m_nu(i)*c**2/(k_B*T_nu0)
do nu_i=1, this%CP%Nu_mass_eigenstates
this%nu_masses(nu_i)= ThermalNuBackground%find_nu_mass_for_rho(this%CP%omnuh2/h2*this%CP%Nu_mass_fractions(nu_i)&
*this%grhocrit/this%grhormass(nu_i))
end do
if (all(this%nu_masses(1:this%CP%Nu_mass_eigenstates)==0)) then
!All density accounted for by massless, so just use massless
this%CP%Omnuh2 = 0
goto 10
end if
!Just prevent divide by zero
this%nu_masses(1:this%CP%Nu_mass_eigenstates) = max(this%nu_masses(1:this%CP%Nu_mass_eigenstates),1e-3_dl)
else
this%nu_masses = 0
end if
call this%CP%DarkEnergy%Init(this)
if (global_error_flag==0) this%tau0=this%TimeOfz(0._dl)
if (global_error_flag==0) then
this%chi0=this%rofChi(this%tau0/this%curvature_radius)
this%scale= this%chi0*this%curvature_radius/this%tau0 !e.g. change l sampling depending on approx peak spacing
if (this%closed .and. this%tau0/this%curvature_radius >3.14) then
call GlobalError('chi >= pi in closed model not supported',error_unsupported_params)
end if
if (WantReion) call this%CP%Reion%Init(this)
if (this%CP%NonLinear/=NonLinear_None .and. .not. back_only) &
call this%CP%NonLinearModel%Init(this)
end if
end if
if (allocated(this%CP%SourceWindows) .and. .not. back_only) then
if (.not. this%CP%WantScalars) then
this%num_redshiftwindows=0
else
this%num_redshiftwindows = size(this%CP%SourceWindows)
end if
else
this%num_redshiftwindows = 0
this%CP%SourceTerms%limber_windows = .false.
endif
if (this%CP%WantScalars .and. this%CP%WantCls .and. this%num_redshiftwindows>0) then
eta_k = this%CP%Max_eta_k
if (allocated(this%Redshift_W)) deallocate(this%Redshift_W)
allocate(this%Redshift_W(this%num_redshiftwindows))
this%num_extra_redshiftwindows = 0
!$OMP PARALLEL DO DEFAULT(SHARED), SCHEDULE(STATIC), PRIVATE(zpeak, sigma_z, zpeakstart, zpeakend, nu_i, Win)
do nu_i = 1, this%num_redshiftwindows
Win => this%Redshift_w(nu_i)
Win%Window => this%CP%SourceWindows(nu_i)%Window
Win%kind = Win%Window%source_type
call Win%Window%GetScales(zpeak, sigma_z, zpeakstart, zpeakend)
if (FeedbackLevel > 1) then
write(*,*) FormatString('Window scales: %d peak: %f, sigma: %f, start:%f, end %f', &
nu_i, zpeak, sigma_z, zpeakstart, zpeakend)
end if
Win%Redshift = zpeak
Win%tau = this%TimeOfz(zpeak, tol=1e-4_dl)
Win%sigma_tau = sigma_z*dtauda(this,1/(1+zpeak))/(1+zpeak)**2
Win%tau_peakstart=this%TimeOfZ(zpeakstart, tol=1e-4_dl)
Win%tau_peakend = this%TimeOfZ(max(0._dl,zpeakend), tol=1e-4_dl)
Win%chi0 = this%tau0-Win%tau
Win%chimin = min(Win%chi0,this%tau0 - this%TimeOfz(max(0.05_dl,zpeakend), tol=1e-4_dl))
!$OMP CRITICAL
this%CP%Max_eta_k = max(this%CP%Max_eta_k, this%tau0*WindowKmaxForL(Win,this%CP,this%CP%max_l))
if (Win%Window%source_type==window_21cm) this%CP%Do21cm = .true.
if (Win%Window%source_type==window_counts .and. P%SourceTerms%counts_lensing) then
this%num_extra_redshiftwindows = this%num_extra_redshiftwindows + 1
Win%mag_index = this%num_extra_redshiftwindows
end if
!$OMP END CRITICAL
end do
if (eta_k /= this%CP%Max_eta_k .and. FeedbackLevel>0) &
write (*,*) 'source max_eta_k: ', this%CP%Max_eta_k,'kmax = ', this%CP%Max_eta_k/this%tau0
end if
if ((this%CP%NonLinear==NonLinear_Lens .or. this%CP%NonLinear==NonLinear_both) .and. &
this%CP%Max_eta_k/this%tau0 > this%CP%Transfer%kmax) then
this%CP%Transfer%kmax =this%CP%Max_eta_k/this%tau0
if (FeedbackLevel > 0) write (*,*) 'kmax changed to ', this%CP%Transfer%kmax
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. calling_again) then
write(*,'("Om_b h^2 = ",f9.6)') P%ombh2
write(*,'("Om_c h^2 = ",f9.6)') P%omch2
write(*,'("Om_nu h^2 = ",f9.6)') P%omnuh2
write(*,'("Om_darkenergy = ",f9.6)') this%Omega_de
write(*,'("Om_K = ",f9.6)') P%omk
write(*,'("Om_m (inc Om_u) = ",f9.6)') (P%ombh2+P%omch2+P%omnuh2)/h2
write(*,'("100 theta (CosmoMC) = ",f9.6)') 100*this%CosmomcTheta()
if (this%CP%Num_Nu_Massive > 0) then
write(*,'("N_eff (total) = ",f9.6)') nu_massless_degeneracy + &
sum(this%CP%Nu_mass_degeneracies(1:this%CP%Nu_mass_eigenstates))
do nu_i=1, this%CP%Nu_mass_eigenstates
conv = k_B*(8*this%grhor/this%grhog/7)**0.25*this%CP%tcmb/eV * &
(this%CP%nu_mass_degeneracies(nu_i)/this%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)")') &
this%CP%nu_mass_numbers(nu_i), this%CP%nu_mass_degeneracies(nu_i), &
this%nu_masses(nu_i),conv*this%nu_masses(nu_i)
end do
end if
end if
end subroutine CAMBdata_SetParams
subroutine CAMBdata_Free(this)
class(CAMBdata) :: this
call Free_ClTransfer(this%CLdata%CTransScal)
call Free_ClTransfer(this%ClData%CTransVec)
call Free_ClTransfer(this%ClData%CTransTens)
call this%MT%Free()
if (allocated(this%CAMB_Pk)) deallocate(this%CAMB_PK)
end subroutine CAMBdata_Free
function CAMBdata_DeltaTime(this, a1,a2, in_tol)
class(CAMBdata) :: this
real(dl) CAMBdata_DeltaTime, atol
real(dl), intent(IN) :: a1,a2
real(dl), optional, intent(in) :: in_tol
atol = PresentDefault(tol/1000/exp(this%CP%Accuracy%AccuracyBoost*this%CP%Accuracy%IntTolBoost-1), in_tol)
CAMBdata_DeltaTime = Integrate_Romberg(this, dtauda,a1,a2,atol)
end function CAMBdata_DeltaTime
subroutine CAMBdata_DeltaTimeArr(this, arr, a1, a2, n, tol)
class(CAMBdata) :: this
integer, intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: a1(n), a2(n)
real(dl), intent(in), optional :: tol
integer i
!$OMP PARALLEL DO DEFAULT(SHARED),SCHEDULE(STATIC)
do i = 1, n
arr(i) = this%DeltaTime(a1(i), a2(i), tol)
end do
end subroutine CAMBdata_DeltaTimeArr
function CAMBdata_TimeOfz(this, z, tol)
class(CAMBdata) :: this
real(dl) CAMBdata_TimeOfz
real(dl), intent(in), optional :: tol
real(dl), intent(IN) :: z
CAMBdata_TimeOfz= this%DeltaTime(0._dl,1._dl/(z+1._dl), tol)
end function CAMBdata_TimeOfz
subroutine CAMBdata_TimeOfzArr(this, arr, z, n, tol)
!z array must be monotonically *decreasing* so times increasing
class(CAMBdata) :: this
integer, intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: z(n)
real(dl), intent(in), optional :: tol
integer i
!$OMP PARALLEL DO DEFAULT(SHARED),SCHEDULE(STATIC)
do i = 1, n
if (i==1) then
arr(i) = this%DeltaTime(0._dl, 1/(1+z(1)), tol)
else
if (z(i) < z(i-1)) then
arr(i) = this%DeltaTime(1/(1+z(i-1)),1/(1+z(i)),tol)
elseif (z(i) < z(i-1) + 1e-6_dl) then
arr(i)=0
else
error stop 'CAMBdata_TimeOfzArr redshifts must be decreasing'
end if
end if
end do
!$OMP END PARALLEL DO
do i = 2, n
arr(i) = arr(i) + arr(i-1)
end do
end subroutine CAMBdata_TimeOfzArr
function CAMBdata_DeltaPhysicalTimeGyr(this, a1,a2, in_tol)
use constants
class(CAMBdata) :: this
real(dl), intent(in) :: a1, a2
real(dl), optional, intent(in) :: in_tol
real(dl) CAMBdata_DeltaPhysicalTimeGyr, atol
atol = PresentDefault(1d-4/exp(this%CP%Accuracy%AccuracyBoost-1), in_tol)
CAMBdata_DeltaPhysicalTimeGyr = Integrate_Romberg(this, dtda,a1,a2,atol)*Mpc/c/Gyr
end function CAMBdata_DeltaPhysicalTimeGyr
subroutine CAMBdata_DeltaPhysicalTimeGyrArr(this, arr, a1, a2, n, tol)
class(CAMBdata) :: this
integer, intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: a1(n), a2(n)
real(dl), intent(in), optional :: tol
integer i
!$OMP PARALLEL DO DEFAULT(SHARED),SCHEDULE(STATIC)
do i = 1, n
arr(i) = this%DeltaPhysicalTimeGyr(a1(i), a2(i), tol)
end do
end subroutine CAMBdata_DeltaPhysicalTimeGyrArr
function CAMBdata_AngularDiameterDistance(this,z)
class(CAMBdata) :: this
!This is the physical (non-comoving) angular diameter distance in Mpc
real(dl) CAMBdata_AngularDiameterDistance
real(dl), intent(in) :: z
CAMBdata_AngularDiameterDistance = this%curvature_radius/(1+z)* &
this%rofchi(this%ComovingRadialDistance(z) /this%curvature_radius)
end function CAMBdata_AngularDiameterDistance
subroutine CAMBdata_AngularDiameterDistanceArr(this, arr, z, n)
class(CAMBdata) :: this
!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 this%ComovingRadialDistanceArr(arr, z, n, 1e-4_dl)
if (this%flat) then
arr = arr/(1+z)
else
do i=1, n
arr(i) = this%curvature_radius/(1+z(i))*this%rofchi(arr(i)/this%curvature_radius)
end do
end if
end subroutine CAMBdata_AngularDiameterDistanceArr
function CAMBdata_AngularDiameterDistance2(this,z1, z2)
! z1 < z2, otherwise returns zero
!From http://www.slac.stanford.edu/~amantz/work/fgas14/#cosmomc
class(CAMBdata) :: this
real(dl) CAMBdata_AngularDiameterDistance2
real(dl), intent(in) :: z1, z2
if (z2 < z1 + 1e-4) then
CAMBdata_AngularDiameterDistance2=0
else
CAMBdata_AngularDiameterDistance2 = this%curvature_radius/(1+z2)* &
this%rofchi( this%DeltaTime(1/(1+z2),1/(1+z1))/this%curvature_radius)
end if
end function CAMBdata_AngularDiameterDistance2
subroutine CAMBdata_AngularDiameterDistance2Arr(this, arr, z1, z2, n)
class(CAMBdata) :: this
integer, intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: z1(n), z2(n)
integer i
!$OMP PARALLEL DO DEFAULT(SHARED),SCHEDULE(STATIC)
do i = 1, n
arr(i) = this%AngularDiameterDistance2(z1(i),z2(i))
end do
!$OMP END PARALLEL DO
end subroutine CAMBdata_AngularDiameterDistance2Arr
function CAMBdata_LuminosityDistance(this,z)
class(CAMBdata) :: this
real(dl) CAMBdata_LuminosityDistance
real(dl), intent(in) :: z
CAMBdata_LuminosityDistance = this%AngularDiameterDistance(z)*(1+z)**2
end function CAMBdata_LuminosityDistance
function CAMBdata_ComovingRadialDistance(this, z)
class(CAMBdata) :: this
real(dl) CAMBdata_ComovingRadialDistance
real(dl), intent(in) :: z
CAMBdata_ComovingRadialDistance = this%DeltaTime(1/(1+z),1._dl)
end function CAMBdata_ComovingRadialDistance
subroutine CAMBdata_ComovingRadialDistanceArr(this, arr, z, n, tol)
!z array must be monotonically increasing
class(CAMBdata) :: this
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) = this%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) = this%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 CAMBdata_ComovingRadialDistanceArr
function CAMBdata_Hofz(this,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
class(CAMBdata) :: this
real(dl) CAMBdata_Hofz, a
real(dl), intent(in) :: z
a = 1/(1+z)
CAMBdata_Hofz = 1/(a**2*dtauda(this,a))
end function CAMBdata_Hofz
subroutine CAMBdata_HofzArr(this, 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
class(CAMBdata) :: this
integer,intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: z(n)
integer i
real(dl) :: a
do i=1, n
a = 1/(1+z(i))
arr(i) = 1/(a**2*dtauda(this,a))
end do
end subroutine CAMBdata_HofzArr
real(dl) function CAMBdata_sound_horizon(this, z)
class(CAMBdata) :: this
real(dl), intent(in) :: z
CAMBdata_sound_horizon = Integrate_Romberg(this,dsound_da_exact,1d-9,1/(z+1),1e-6_dl)
end function CAMBdata_sound_horizon
subroutine CAMBdata_sound_horizon_zArr(this,arr, z,n)
class(CAMBdata) :: this
integer,intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: z(n)
integer i
!$OMP PARALLEL DO DEFAULT(SHARED), SCHEDULE(STATIC), IF(n>4)
do i=1,n
arr(i) = this%sound_horizon(z(i))
end do
end subroutine CAMBdata_sound_horizon_zArr
subroutine CAMBdata_RedshiftAtTimeArr(this, arr, tau, n)
class(CAMBdata) :: this
integer,intent(in) :: n
real(dl), intent(out) :: arr(n)
real(dl), intent(in) :: tau(n)
integer i
real(dl) om
if (this%ThermoData%ScaleFactorAtTime%n==0) &
call GlobalError('RedshiftAtTimeArr: background history not calculated', error_unsupported_params)
if (global_error_flag/=0) return
!$OMP PARALLEL DO DEFAULT(SHARED), private(om, i)
do i=1, n
if (tau(i) < this%ThermoData%tauminn*1.1) then
om = (this%grhob+this%grhoc)/&
sqrt(3*(this%grhog+sum(this%grhormass(1:this%CP%Nu_mass_eigenstates))+this%grhornomass))
arr(i) = 1/(this%adotrad*tau(i)*(1+om*tau(i)/4))-1
else
arr(i) = 1/this%ThermoData%ScaleFactorAtTime%Value(tau(i))-1
end if
end do
end subroutine CAMBdata_RedshiftAtTimeArr
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 CAMBdata_BAO_D_v(this,z)
class(CAMBdata) :: this
real(dl), intent(IN) :: z
CAMBdata_BAO_D_v = BAO_D_v_from_DA_H(z,this%AngularDiameterDistance(z), this%Hofz(z))
end function CAMBdata_BAO_D_v
function CAMBdata_get_zstar(this)
class(CAMBdata) :: this
real(dl) CAMBdata_get_zstar
real(dl) z_scale
call this%CP%Recomb%Init(this)
z_scale = COBE_CMBTemp/this%CP%TCMB
CAMBdata_get_zstar=this%binary_search(noreion_optdepth, 1.d0, 700.d0*z_scale, &
2000.d0*z_scale, 1d-3*z_scale,100.d0*z_scale,5000.d0*z_scale)
end function CAMBdata_get_zstar
function CAMBdata_CosmomcTheta(this)
class(CAMBdata) :: this
real(dl) zstar, astar, atol, rs, DA
real(dl) CAMBdata_CosmomcTheta
real(dl) ombh2, omdmh2
ombh2 = this%CP%ombh2
omdmh2 = (this%CP%omch2+this%CP%omnuh2)
!!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 = Integrate_Romberg(this,dsound_da_approx,1d-8,astar,atol)
DA = this%AngularDiameterDistance(zstar)/astar
CAMBdata_CosmomcTheta = rs/DA
end function CAMBdata_CosmomcTheta
subroutine CAMBdata_GetBackgroundDensities(this, n, a_arr, densities)
! return array of 8*pi*G*rho*a**4 for each species
class(CAMBdata) :: this
integer, intent(in) :: n
real(dl), intent(in) :: a_arr(n)
real(dl) :: grhov_t, rhonu, grhonu, a
real(dl), intent(out) :: densities(8,n)
integer nu_i,i
do i=1, n
a = a_arr(i)
call this%CP%DarkEnergy%BackgroundDensityAndPressure(this%grhov, a, grhov_t)
grhonu = 0
if (this%CP%Num_Nu_massive /= 0) then
!Get massive neutrino density relative to massless
do nu_i = 1, this%CP%nu_mass_eigenstates
call ThermalNuBackground%rho(a * this%nu_masses(nu_i), rhonu)
grhonu = grhonu + rhonu * this%grhormass(nu_i)
end do
end if
densities(2,i) = this%grhok * a**2
densities(3,i) = this%grhoc * a
densities(4,i) = this%grhob * a
densities(5,i) = this%grhog
densities(6,i) = this%grhornomass
densities(7,i) = grhonu
densities(8,i) = grhov_t*a**2
densities(1,i) = sum(densities(2:8,i))
end do
end subroutine CAMBdata_GetBackgroundDensities
integer function CAMBdata_get_lmax_lensed(this)
class(CAMBdata) :: this
CAMBdata_get_lmax_lensed = this%CLdata%lmax_lensed
end function CAMBdata_get_lmax_lensed
!JD 08/13 New function for nonlinear lensing of CMB + MPK compatibility
!Build master redshift array from array of desired Nonlinear lensing (NLL)
!redshifts and an array of desired Power spectrum (PK) redshifts.
!At the same time fill arrays for NLL and PK that indicate indices
!of their desired redshifts in the master redshift array.
!Finally define number of redshifts in master array. This is usually given by:
!P%num_redshifts = P%PK_num_redshifts + NLL_num_redshifts - 1. The -1 comes
!from the fact that z=0 is in both arrays (when non-linear is on)
subroutine GetComputedPKRedshifts(this, Params,eta_k_max)
use MpiUtils, only : MpiStop
class(CAMBdata) :: this
Type(CAMBParams) :: Params
integer i, iPK, iNLL
real(dl), parameter :: tol = 1.d-5
real(dl) maxRedshift, NL_Boost
integer :: NLL_num_redshifts
real(dl), allocatable :: NLL_redshifts(:), redshifts(:)
!Sources, but unused currently
real(dl), intent(in), optional :: eta_k_max
NLL_num_redshifts = 0
this%needs_good_pk_sampling = .false.
associate(P => Params%Transfer)
if ((Params%NonLinear==NonLinear_lens .or. Params%NonLinear==NonLinear_both) .and. &
(Params%DoLensing .or. this%num_redshiftwindows > 0)) then
! Want non-linear lensing or other sources
this%needs_good_pk_sampling = .false.
NL_Boost = Params%Accuracy%AccuracyBoost*Params%Accuracy%NonlinSourceBoost
if (Params%Do21cm) then
!Sources
if (maxval(this%Redshift_w(1:this%num_redshiftwindows)%Redshift) &
/= minval(this%Redshift_w(1:this%num_redshiftwindows)%Redshift)) &
stop 'Non-linear 21cm currently only for narrow window at one redshift'