CosmicFish  Reference documentation for version 1.0
Looking into future Cosmology
006_cosmicfish_camb_interface.f90
Go to the documentation of this file.
1 !----------------------------------------------------------------------------------------
2 !
3 ! This file is part of CosmicFish.
4 !
5 ! Copyright (C) 2015-2016 by the CosmicFish authors
6 !
7 ! The CosmicFish code is free software;
8 ! You can use it, redistribute it, and/or modify it under the terms
9 ! of the GNU General Public License as published by the Free Software Foundation;
10 ! either version 3 of the License, or (at your option) any later version.
11 ! The full text of the license can be found in the file LICENSE at
12 ! the top level of the CosmicFish distribution.
13 !
14 !----------------------------------------------------------------------------------------
15 
19 
20 !----------------------------------------------------------------------------------------
22 
24 
26 
27  use precision
29  use modelparams
30  use camb
31  use lambdageneral
32 
33  implicit none
34 
35  private
36 
39 
40 contains
41 
42  ! ---------------------------------------------------------------------------------------------
45  subroutine num_param_fisher(P, FP, num_param)
46 
47  implicit none
48 
49  Type(cambparams) , intent(in) :: P
50  Type(cosmicfish_params), intent(in) :: FP
51  integer , intent(out) :: num_param
52 
53  integer :: i
54 
55  num_param = 0
56 
57  ! standard parameters:
58  if ( fp%fisher_par%want_ombh2 ) num_param = num_param + 1 ! ombh2
59  if ( fp%fisher_par%want_omch2 ) num_param = num_param + 1 ! omch2
60  if ( fp%fisher_par%want_omnuh2 ) num_param = num_param + 1 ! omnuh2
61  if ( fp%fisher_par%want_hubble ) num_param = num_param + 1 ! small hubble
62  if ( fp%fisher_par%want_helium_fraction ) num_param = num_param + 1 ! helium abundance
63  if ( fp%fisher_par%want_massless ) num_param = num_param + 1 ! massless neutrinos
64  ! primordial spectrum parameters:
65  if ( fp%fisher_par%want_scalar_amp ) num_param = num_param + 1 ! scalar_amp
66  if ( fp%fisher_par%want_scalar_spectral_index ) num_param = num_param + 1 ! scalar_spectral_index
67  if ( fp%fisher_par%want_scalar_nrun ) num_param = num_param + 1 ! scalar_nrun
68  if ( fp%fisher_par%want_tensor_spectral_index ) num_param = num_param + 1 ! tensor_spectral_index
69  if ( fp%fisher_par%want_initial_ratio ) num_param = num_param + 1 ! initial_ratio
70  ! reionization:
71  if ( fp%fisher_par%want_re_optical_depth ) num_param = num_param + 1 ! reionization optical depth
72 
73  ! Nuisance parameters:
74  ! bias in every counts redshift windows:
75  if ( fp%fisher_par%want_bias ) then
76  do i = 1, num_redshiftwindows
77  if ( redshift_w(i)%kind == window_counts ) then
78  num_param = num_param + 1
79  end if
80  end do
81  end if
82  ! nuissance parameters for supernovae
83  if ( fp%fisher_par%want_alpha_SN ) num_param = num_param + 1 ! alpha_SN
84  if ( fp%fisher_par%want_beta_SN ) num_param = num_param + 1 ! beta_SN
85  if ( fp%fisher_par%want_M0_SN ) num_param = num_param + 1 ! M0_SN
86 
87 #ifdef COSMICFISH_CAMB
88  ! ppf parameters
89  if ( fp%fisher_par%want_w0_ppf ) num_param = num_param + 1 ! w0_ppf
90  if ( fp%fisher_par%want_wa_ppf ) num_param = num_param + 1 ! wa_ppf
91  if ( fp%fisher_par%want_cs_ppf ) num_param = num_param + 1 ! cs_ppf
92 #endif
93 
94 #ifdef COSMICFISH_EFTCAMB
95  ! EFTCAMB parameters
96  ! w_DE parameters:
97  if ( ( p%EFTflag == 1 .or. p%EFTflag == 2 .or. p%EFTflag == 3 ) .and. p%EFTwDE /= 0 ) then
98  if ( p%EFTwDE == 1 ) then ! DE with constant Eos
99  num_param = num_param + 1 ! EFTw0
100  else if ( p%EFTwDE == 2 ) then ! CPL parametrization
101  num_param = num_param + 2 ! EFTw0, EFTwa
102  else if ( p%EFTwDE == 3 ) then ! JBP parametrization
103  num_param = num_param + 3 ! EFTw0, EFTwa, EFTwn
104  else if ( p%EFTwDE == 4 ) then ! turning point parametrization
105  num_param = num_param + 3 ! EFTw0, EFTwa, EFTwat
106  else if ( p%EFTwDE == 5 ) then ! Taylor expansion
107  num_param = num_param + 4 ! EFTw0, EFTwa, EFTw2, EFTw3
108  else if ( p%EFTwDE == 6 ) then ! User defined
109  num_param = num_param + 1
110  end if
111  end if
112  ! EFT functions parameters:
113  if ( p%EFTflag == 1 ) then ! Pure EFT models
114  ! EFTOmega:
115  if ( p%PureEFTmodelOmega == 1 ) then
116  num_param = num_param + 1
117  else if ( p%PureEFTmodelOmega == 2 ) then
118  num_param = num_param + 1
119  else if ( p%PureEFTmodelOmega == 3 ) then
120  num_param = num_param + 2
121  else if ( p%PureEFTmodelOmega == 4 ) then
122  num_param = num_param + 2
123  else if ( p%PureEFTmodelOmega == 5 ) then
124  num_param = num_param + 1
125  end if
126  ! EFTGamma1
127  if ( p%PureEFTmodelGamma1 == 1 ) then
128  num_param = num_param + 1
129  else if ( p%PureEFTmodelGamma1 == 2 ) then
130  num_param = num_param + 1
131  else if ( p%PureEFTmodelGamma1 == 3 ) then
132  num_param = num_param + 2
133  else if ( p%PureEFTmodelGamma1 == 4 ) then
134  num_param = num_param + 2
135  else if ( p%PureEFTmodelGamma1 == 5 ) then
136  num_param = num_param + 1
137  end if
138  ! EFTGamma2
139  if ( p%PureEFTmodelGamma2 == 1 ) then
140  num_param = num_param + 1
141  else if ( p%PureEFTmodelGamma2 == 2 ) then
142  num_param = num_param + 1
143  else if ( p%PureEFTmodelGamma2 == 3 ) then
144  num_param = num_param + 2
145  else if ( p%PureEFTmodelGamma2 == 4 ) then
146  num_param = num_param + 2
147  else if ( p%PureEFTmodelGamma2 == 5 ) then
148  num_param = num_param + 1
149  end if
150  ! EFTGamma3
151  if ( p%PureEFTmodelGamma3 == 1 ) then
152  num_param = num_param + 1
153  else if ( p%PureEFTmodelGamma3 == 2 ) then
154  num_param = num_param + 1
155  else if ( p%PureEFTmodelGamma3 == 3 ) then
156  num_param = num_param + 2
157  else if ( p%PureEFTmodelGamma3 == 4 ) then
158  num_param = num_param + 2
159  else if ( p%PureEFTmodelGamma3 == 5 ) then
160  num_param = num_param + 1
161  end if
162  ! EFTGamma4
163  if ( p%PureEFTmodelGamma4 == 1 ) then
164  num_param = num_param + 1
165  else if ( p%PureEFTmodelGamma4 == 2 ) then
166  num_param = num_param + 1
167  else if ( p%PureEFTmodelGamma4 == 3 ) then
168  num_param = num_param + 2
169  else if ( p%PureEFTmodelGamma4 == 4 ) then
170  num_param = num_param + 2
171  else if ( p%PureEFTmodelGamma4 == 5 ) then
172  num_param = num_param + 1
173  end if
174  ! EFTGamma5
175  if ( p%PureEFTmodelGamma5 == 1 ) then
176  num_param = num_param + 1
177  else if ( p%PureEFTmodelGamma5 == 2 ) then
178  num_param = num_param + 1
179  else if ( p%PureEFTmodelGamma5 == 3 ) then
180  num_param = num_param + 2
181  else if ( p%PureEFTmodelGamma5 == 4 ) then
182  num_param = num_param + 2
183  else if ( p%PureEFTmodelGamma5 == 5 ) then
184  num_param = num_param + 1
185  end if
186  ! EFTGamma6
187  if ( p%PureEFTmodelGamma6 == 1 ) then
188  num_param = num_param + 1
189  else if ( p%PureEFTmodelGamma6 == 2 ) then
190  num_param = num_param + 1
191  else if ( p%PureEFTmodelGamma6 == 3 ) then
192  num_param = num_param + 2
193  else if ( p%PureEFTmodelGamma6 == 4 ) then
194  num_param = num_param + 2
195  else if ( p%PureEFTmodelGamma6 == 5 ) then
196  num_param = num_param + 1
197  end if
198  else if ( p%EFTflag == 2 ) then ! designer EFT models
199  if ( p%DesignerEFTmodel == 1 ) then
200  num_param = num_param + 1 ! B0
201  else if ( p%DesignerEFTmodel == 2 ) then
202  num_param = num_param + 0 ! mc 5e models
203  end if
204  else if ( p%EFTflag == 3 ) then ! EFT alternative parametrizations
205  if ( p%AltParEFTmodel == 1 ) then
206  if ( p%RPHmassPmodel == 1 ) then
207  num_param = num_param + 1
208  else if ( p%RPHmassPmodel == 2 ) then
209  num_param = num_param + 2
210  end if
211  if ( p%RPHkineticitymodel == 1 ) then
212  num_param = num_param + 1
213  else if ( p%RPHkineticitymodel == 2 ) then
214  num_param = num_param + 2
215  end if
216  if ( p%RPHbraidingmodel == 1 ) then
217  num_param = num_param + 1
218  else if ( p%RPHbraidingmodel == 2 ) then
219  num_param = num_param + 2
220  end if
221  if ( p%RPHtensormodel == 1 ) then
222  num_param = num_param + 1
223  else if ( p%RPHtensormodel == 2 ) then
224  num_param = num_param + 2
225  end if
226  end if
227  else if ( p%EFTflag == 4 ) then ! Full mapping EFT model
228  if ( p%FullMappingEFTmodel == 1 ) then
229  if ( p%HoravaSolarSystem .eqv. .true. ) then
230  num_param = num_param + 2
231  else
232  num_param = num_param + 3
233  end if
234  end if
235  end if
236 #endif
237 
238 #ifdef COSMICFISH_MGCAMB
239  if ( p%MGC_model == 1 ) then
240  num_param = num_param + 5
241  else if ( p%MGC_model == 2 ) then
242  num_param = num_param + 2
243  else if ( p%MGC_model == 3 ) then
244  num_param = num_param + 3
245  else if ( p%MGC_model == 4 ) then
246  num_param = num_param + 1
247  else if ( p%MGC_model ==5 ) then
248  num_param = num_param + 3
249  else if ( p%MGC_model ==6 ) then
250  num_param = num_param + 1
251  else if ( p%MGC_model == 7 ) then
252  num_param = num_param + 3
253  else if ( p%MGC_model == 8 ) then
254  num_param = num_param + 4
255  else if ( p%MGC_model == 9 ) then
256  num_param = num_param + 2
257  else if ( p%MGC_model ==10 ) then
258  num_param = num_param + 2
259  end if
260 #endif
261 
262  end subroutine num_param_fisher
263 
264  ! ---------------------------------------------------------------------------------------------
266  subroutine fisher_param_names(P, FP, param_number, param_name, param_name_latex)
268  implicit none
269 
270  Type(cambparams) , intent(in) :: P
271  Type(cosmicfish_params), intent(in) :: FP
272  integer , intent(in) :: param_number
273  character(*) , intent(out) :: param_name
274  character(*) , intent(out), optional :: param_name_latex
275 
276  integer :: i, num_param
277  character(len=20) str
278 
279  num_param = 0
280 
281  ! standard parameters:
282  if ( fp%fisher_par%want_ombh2 ) num_param = num_param + 1 ! ombh2
283  if ( num_param == param_number ) then; param_name = 'omegabh2'; if ( present(param_name_latex) ) param_name_latex = '\Omega_b h^2'; return ;
284  end if
285 
286  if ( fp%fisher_par%want_omch2 ) num_param = num_param + 1 ! omch2
287  if ( num_param == param_number ) then; param_name = 'omegach2'; if ( present(param_name_latex) ) param_name_latex = '\Omega_c h^2'; return ;
288  end if
289 
290  if ( fp%fisher_par%want_omnuh2 ) num_param = num_param + 1 ! omnuh2
291  if ( num_param == param_number ) then; param_name = 'omeganuh2'; if ( present(param_name_latex) ) param_name_latex = '\Omega_\nu h^2'; return ;
292  end if
293 
294  if ( fp%fisher_par%want_hubble ) num_param = num_param + 1 ! omegav or small hubble
295  if ( num_param == param_number ) then; param_name = 'h'; if ( present(param_name_latex) ) param_name_latex = 'h'; return ;
296  end if
297 
298  if ( fp%fisher_par%want_helium_fraction ) num_param = num_param + 1 ! helium abundance
299  if ( num_param == param_number ) then; param_name = 'yhe'; if ( present(param_name_latex) ) param_name_latex = 'Y_{He}'; return ;
300  end if
301 
302  if ( fp%fisher_par%want_massless ) num_param = num_param + 1 ! massless neutrinos
303  if ( num_param == param_number ) then; param_name = 'numless'; if ( present(param_name_latex) ) param_name_latex = '\nu_{mless}'; return ;
304  end if
305 
306  ! primordial spectrum parameters:
307  if ( fp%fisher_par%want_scalar_amp ) num_param = num_param + 1 ! scalar_amp
308  if ( num_param == param_number ) then; param_name = 'logA'; if ( present(param_name_latex) ) param_name_latex = '{\rm{ln}}(10^{10} A_s)'; return ;
309  end if
310 
311  if ( fp%fisher_par%want_scalar_spectral_index ) num_param = num_param + 1 ! scalar_spectral_index
312  if ( num_param == param_number ) then; param_name = 'ns'; if ( present(param_name_latex) ) param_name_latex = 'n_s'; return ;
313  end if
314 
315  if ( fp%fisher_par%want_scalar_nrun ) num_param = num_param + 1 ! scalar_nrun
316  if ( num_param == param_number ) then; param_name = 'nrun'; if ( present(param_name_latex) ) param_name_latex = 'n_{\rm run}'; return ;
317  end if
318 
319  if ( fp%fisher_par%want_tensor_spectral_index ) num_param = num_param + 1 ! tensor_spectral_index
320  if ( num_param == param_number ) then; param_name = 'nt'; if ( present(param_name_latex) ) param_name_latex = 'n_t'; return ;
321  end if
322 
323  if ( fp%fisher_par%want_initial_ratio ) num_param = num_param + 1 ! initial_ratio
324  if ( num_param == param_number ) then; param_name = 'r'; if ( present(param_name_latex) ) param_name_latex = 'r'; return ;
325  end if
326 
327  ! reionization:
328  if ( fp%fisher_par%want_re_optical_depth ) num_param = num_param + 1 ! reionization optical depth
329  if ( num_param == param_number ) then; param_name = 'tau'; if ( present(param_name_latex) ) param_name_latex = '\tau'; return ;
330  end if
331 
332  ! Nuisance parameters:
333  ! bias in every counts redshift windows:
334  if ( fp%fisher_par%want_bias ) then
335  do i = 1, num_redshiftwindows
336  if ( redshift_w(i)%kind == window_counts ) then
337  num_param = num_param + 1
338  write (str, *) i
339  str = adjustl(str)
340  if ( num_param == param_number ) then; param_name = 'Bias_W_'//trim(str) ; if ( present(param_name_latex) ) param_name_latex = 'b_{'//trim(str)//'}'; return ;
341  end if
342  end if
343  end do
344  end if
345 
346  ! nuissance parameters for supernovae
347  if ( fp%fisher_par%want_alpha_SN ) num_param = num_param + 1 ! alpha_SN
348  if ( num_param == param_number ) then; param_name = 'alpha_SN'; if ( present(param_name_latex) ) param_name_latex = '\alpha_{\rm SN}'; return ;
349  end if
350  if ( fp%fisher_par%want_beta_SN ) num_param = num_param + 1 ! beta_SN
351  if ( num_param == param_number ) then; param_name = 'beta_SN'; if ( present(param_name_latex) ) param_name_latex = '\beta_{\rm SN}'; return ;
352  end if
353  if ( fp%fisher_par%want_M0_SN ) num_param = num_param + 1 ! M0_SN
354  if ( num_param == param_number ) then; param_name = 'M0_SN'; if ( present(param_name_latex) ) param_name_latex = 'M_0^{\rm SN}'; return ;
355  end if
356 
357 #ifdef COSMICFISH_CAMB
358  ! ppf parameters
359  if ( fp%fisher_par%want_w0_ppf ) num_param = num_param + 1 ! w0_ppf
360  if ( num_param == param_number ) then; param_name = 'w0_ppf'; if ( present(param_name_latex) ) param_name_latex = 'w_{0}^{\rm ppf}'; return ;
361  end if
362  if ( fp%fisher_par%want_wa_ppf ) num_param = num_param + 1 ! wa_ppf
363  if ( num_param == param_number ) then; param_name = 'wa_ppf'; if ( present(param_name_latex) ) param_name_latex = 'w_{a}^{\rm ppf}'; return ;
364  end if
365  if ( fp%fisher_par%want_cs_ppf ) num_param = num_param + 1 ! cs_ppf
366  if ( num_param == param_number ) then; param_name = 'cs_ppf'; if ( present(param_name_latex) ) param_name_latex = 'c_{s}^{\rm ppf}'; return ;
367  end if
368 #endif
369 
370 #ifdef COSMICFISH_EFTCAMB
371  ! EFTCAMB parameters
372  ! w_DE parameters:
373  if ( ( p%EFTflag == 1 .or. p%EFTflag == 2 .or. p%EFTflag == 3 ) .and. p%EFTwDE /= 0 ) then
374  if ( p%EFTwDE == 1 ) then ! DE with constant Eos
375  num_param = num_param + 1 ! EFTw0
376  if ( num_param == param_number ) then; param_name = 'EFTw0' ; if ( present(param_name_latex) ) param_name_latex = 'w_0'; return ;
377  end if
378  else if ( p%EFTwDE == 2 ) then ! CPL parametrization
379  num_param = num_param + 1 ! EFTw0
380  if ( num_param == param_number ) then; param_name = 'EFTw0' ; if ( present(param_name_latex) ) param_name_latex = 'w_0'; return ;
381  end if
382  num_param = num_param + 1 ! EFTwa
383  if ( num_param == param_number ) then; param_name = 'EFTwa' ; if ( present(param_name_latex) ) param_name_latex = 'w_a'; return ;
384  end if
385  else if ( p%EFTwDE == 3 ) then ! JBP parametrization
386  num_param = num_param + 1 ! EFTw0
387  if ( num_param == param_number ) then; param_name = 'EFTw0' ; if ( present(param_name_latex) ) param_name_latex = 'w_0'; return ;
388  end if
389  num_param = num_param + 1 ! EFTwa
390  if ( num_param == param_number ) then; param_name = 'EFTwa' ; if ( present(param_name_latex) ) param_name_latex = 'w_a'; return ;
391  end if
392  num_param = num_param + 1 ! EFTwn
393  if ( num_param == param_number ) then; param_name = 'EFTwn' ; if ( present(param_name_latex) ) param_name_latex = 'w_n'; return ;
394  end if
395  else if ( p%EFTwDE == 4 ) then ! turning point parametrization
396  num_param = num_param + 1 ! EFTw0
397  if ( num_param == param_number ) then; param_name = 'EFTw0' ; if ( present(param_name_latex) ) param_name_latex = 'w_0'; return ;
398  end if
399  num_param = num_param + 1 ! EFTwa
400  if ( num_param == param_number ) then; param_name = 'EFTwa' ; if ( present(param_name_latex) ) param_name_latex = 'w_a'; return ;
401  end if
402  num_param = num_param + 1 ! EFTwat
403  if ( num_param == param_number ) then; param_name = 'EFTwat' ; if ( present(param_name_latex) ) param_name_latex = 'wa_t'; return ;
404  end if
405  else if ( p%EFTwDE == 5 ) then ! Taylor expansion
406  num_param = num_param + 1 ! EFTw0
407  if ( num_param == param_number ) then; param_name = 'EFTw0' ; if ( present(param_name_latex) ) param_name_latex = 'w_0'; return ;
408  end if
409  num_param = num_param + 1 ! EFTwa
410  if ( num_param == param_number ) then; param_name = 'EFTwa' ; if ( present(param_name_latex) ) param_name_latex = 'w_a'; return ;
411  end if
412  num_param = num_param + 1 ! EFTw2
413  if ( num_param == param_number ) then; param_name = 'EFTw2' ; if ( present(param_name_latex) ) param_name_latex = 'w_2'; return ;
414  end if
415  num_param = num_param + 1 ! EFTw3
416  if ( num_param == param_number ) then; param_name = 'EFTw3' ; if ( present(param_name_latex) ) param_name_latex = 'w_3'; return ;
417  end if
418  else if ( p%EFTwDE == 6 ) then ! User defined
419 
420  end if
421  end if
422 
423  ! EFT functions parameters:
424  if ( p%EFTflag == 1 ) then ! Pure EFT models
425  ! EFTOmega:
426  if ( p%PureEFTmodelOmega == 1 ) then
427  num_param = num_param + 1
428  if ( num_param == param_number ) then; param_name = 'EFTOmega0' ; if ( present(param_name_latex) ) param_name_latex = '\Omega_0^{\rm EFT}'; return ;
429  end if
430  else if ( p%PureEFTmodelOmega == 2 ) then
431  num_param = num_param + 1
432  if ( num_param == param_number ) then; param_name = 'EFTOmega0' ; if ( present(param_name_latex) ) param_name_latex = '\Omega_0^{\rm EFT}'; return ;
433  end if
434  else if ( p%PureEFTmodelOmega == 3 ) then
435  num_param = num_param + 1
436  if ( num_param == param_number ) then; param_name = 'EFTOmega0' ; if ( present(param_name_latex) ) param_name_latex = '\Omega_0^{\rm EFT}'; return ;
437  end if
438  num_param = num_param + 1
439  if ( num_param == param_number ) then; param_name = 'EFTOmegaExp' ; if ( present(param_name_latex) ) param_name_latex = 'n^{\rm EFT}'; return ;
440  end if
441  else if ( p%PureEFTmodelOmega == 4 ) then
442  num_param = num_param + 1
443  if ( num_param == param_number ) then; param_name = 'EFTOmega0' ; if ( present(param_name_latex) ) param_name_latex = '\Omega_0^{\rm EFT}'; return ;
444  end if
445  num_param = num_param + 1
446  if ( num_param == param_number ) then; param_name = 'EFTOmegaExp' ; if ( present(param_name_latex) ) param_name_latex = 'n^{\rm EFT}'; return ;
447  end if
448  else if ( p%PureEFTmodelOmega == 5 ) then
449 
450  end if
451  ! EFTGamma1
452  if ( p%PureEFTmodelGamma1 == 1 ) then
453  num_param = num_param + 1
454  if ( num_param == param_number ) then; param_name = 'EFTGamma10' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(1) {\rm EFT}}'; return ;
455  end if
456  else if ( p%PureEFTmodelGamma1 == 2 ) then
457  num_param = num_param + 1
458  if ( num_param == param_number ) then; param_name = 'EFTGamma10' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(1) {\rm EFT}}'; return ;
459  end if
460  else if ( p%PureEFTmodelGamma1 == 3 ) then
461  num_param = num_param + 1
462  if ( num_param == param_number ) then; param_name = 'EFTGamma10' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(1) {\rm EFT}}'; return ;
463  end if
464  num_param = num_param + 1
465  if ( num_param == param_number ) then; param_name = 'EFTGamma1Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(1) {\rm EFT}}'; return ;
466  end if
467  else if ( p%PureEFTmodelGamma1 == 4 ) then
468  num_param = num_param + 1
469  if ( num_param == param_number ) then; param_name = 'EFTGamma10' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(1) {\rm EFT}}'; return ;
470  end if
471  num_param = num_param + 1
472  if ( num_param == param_number ) then; param_name = 'EFTGamma1Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(1) {\rm EFT}}'; return ;
473  end if
474  else if ( p%PureEFTmodelGamma1 == 5 ) then
475 
476  end if
477  ! EFTGamma2
478  if ( p%PureEFTmodelGamma2 == 1 ) then
479  num_param = num_param + 1
480  if ( num_param == param_number ) then; param_name = 'EFTGamma20' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(2) {\rm EFT}}'; return ;
481  end if
482  else if ( p%PureEFTmodelGamma2 == 2 ) then
483  num_param = num_param + 1
484  if ( num_param == param_number ) then; param_name = 'EFTGamma20' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(2) {\rm EFT}}'; return ;
485  end if
486  else if ( p%PureEFTmodelGamma2 == 3 ) then
487  num_param = num_param + 1
488  if ( num_param == param_number ) then; param_name = 'EFTGamma20' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(2) {\rm EFT}}'; return ;
489  end if
490  num_param = num_param + 1
491  if ( num_param == param_number ) then; param_name = 'EFTGamma2Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(2) {\rm EFT}}'; return ;
492  end if
493  else if ( p%PureEFTmodelGamma2 == 4 ) then
494  num_param = num_param + 1
495  if ( num_param == param_number ) then; param_name = 'EFTGamma20' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(2) {\rm EFT}}'; return ;
496  end if
497  num_param = num_param + 1
498  if ( num_param == param_number ) then; param_name = 'EFTGamma2Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(2) {\rm EFT}}'; return ;
499  end if
500  else if ( p%PureEFTmodelGamma2 == 5 ) then
501 
502  end if
503  ! EFTGamma3
504  if ( p%PureEFTmodelGamma3 == 1 ) then
505  num_param = num_param + 1
506  if ( num_param == param_number ) then; param_name = 'EFTGamma30' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(3) {\rm EFT}}'; return ;
507  end if
508  else if ( p%PureEFTmodelGamma3 == 2 ) then
509  num_param = num_param + 1
510  if ( num_param == param_number ) then; param_name = 'EFTGamma30' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(3) {\rm EFT}}'; return ;
511  end if
512  else if ( p%PureEFTmodelGamma3 == 3 ) then
513  num_param = num_param + 1
514  if ( num_param == param_number ) then; param_name = 'EFTGamma30' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(3) {\rm EFT}}'; return ;
515  end if
516  num_param = num_param + 1
517  if ( num_param == param_number ) then; param_name = 'EFTGamma3Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(3) {\rm EFT}}'; return ;
518  end if
519  else if ( p%PureEFTmodelGamma3 == 4 ) then
520  num_param = num_param + 1
521  if ( num_param == param_number ) then; param_name = 'EFTGamma30' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(3) {\rm EFT}}'; return ;
522  end if
523  num_param = num_param + 1
524  if ( num_param == param_number ) then; param_name = 'EFTGamma3Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(3) {\rm EFT}}'; return ;
525  end if
526  else if ( p%PureEFTmodelGamma3 == 5 ) then
527 
528  end if
529  ! EFTGamma4
530  if ( p%PureEFTmodelGamma4 == 1 ) then
531  num_param = num_param + 1
532  if ( num_param == param_number ) then; param_name = 'EFTGamma40' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(4) {\rm EFT}}'; return ;
533  end if
534  else if ( p%PureEFTmodelGamma4 == 2 ) then
535  num_param = num_param + 1
536  if ( num_param == param_number ) then; param_name = 'EFTGamma40' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(4) {\rm EFT}}'; return ;
537  end if
538  else if ( p%PureEFTmodelGamma4 == 3 ) then
539  num_param = num_param + 1
540  if ( num_param == param_number ) then; param_name = 'EFTGamma40' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(4) {\rm EFT}}'; return ;
541  end if
542  num_param = num_param + 1
543  if ( num_param == param_number ) then; param_name = 'EFTGamma4Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(4) {\rm EFT}}'; return ;
544  end if
545  else if ( p%PureEFTmodelGamma4 == 4 ) then
546  num_param = num_param + 1
547  if ( num_param == param_number ) then; param_name = 'EFTGamma40' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(4) {\rm EFT}}'; return ;
548  end if
549  num_param = num_param + 1
550  if ( num_param == param_number ) then; param_name = 'EFTGamma4Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(4) {\rm EFT}}'; return ;
551  end if
552  else if ( p%PureEFTmodelGamma4 == 5 ) then
553 
554  end if
555  ! EFTGamma5
556  if ( p%PureEFTmodelGamma5 == 1 ) then
557  num_param = num_param + 1
558  if ( num_param == param_number ) then; param_name = 'EFTGamma50' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(5) {\rm EFT}}'; return ;
559  end if
560  else if ( p%PureEFTmodelGamma5 == 2 ) then
561  num_param = num_param + 1
562  if ( num_param == param_number ) then; param_name = 'EFTGamma50' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(5) {\rm EFT}}'; return ;
563  end if
564  else if ( p%PureEFTmodelGamma5 == 3 ) then
565  num_param = num_param + 1
566  if ( num_param == param_number ) then; param_name = 'EFTGamma50' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(5) {\rm EFT}}'; return ;
567  end if
568  num_param = num_param + 1
569  if ( num_param == param_number ) then; param_name = 'EFTGamma5Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(5) {\rm EFT}}'; return ;
570  end if
571  else if ( p%PureEFTmodelGamma5 == 4 ) then
572  num_param = num_param + 1
573  if ( num_param == param_number ) then; param_name = 'EFTGamma50' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(5) {\rm EFT}}'; return ;
574  end if
575  num_param = num_param + 1
576  if ( num_param == param_number ) then; param_name = 'EFTGamma5Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(5) {\rm EFT}}'; return ;
577  end if
578  else if ( p%PureEFTmodelGamma5 == 5 ) then
579 
580  end if
581  ! EFTGamma6
582  if ( p%PureEFTmodelGamma6 == 1 ) then
583  num_param = num_param + 1
584  if ( num_param == param_number ) then; param_name = 'EFTGamma60' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(6) {\rm EFT}}'; return ;
585  end if
586  else if ( p%PureEFTmodelGamma6 == 2 ) then
587  num_param = num_param + 1
588  if ( num_param == param_number ) then; param_name = 'EFTGamma60' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(6) {\rm EFT}}'; return ;
589  end if
590  else if ( p%PureEFTmodelGamma6 == 3 ) then
591  num_param = num_param + 1
592  if ( num_param == param_number ) then; param_name = 'EFTGamma60' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(6) {\rm EFT}}'; return ;
593  end if
594  num_param = num_param + 1
595  if ( num_param == param_number ) then; param_name = 'EFTGamma6Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(6) {\rm EFT}}'; return ;
596  end if
597  else if ( p%PureEFTmodelGamma6 == 4 ) then
598  num_param = num_param + 1
599  if ( num_param == param_number ) then; param_name = 'EFTGamma60' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_0^{(6) {\rm EFT}}'; return ;
600  end if
601  num_param = num_param + 1
602  if ( num_param == param_number ) then; param_name = 'EFTGamma6Exp' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm exp}^{(6) {\rm EFT}}'; return ;
603  end if
604  else if ( p%PureEFTmodelGamma6 == 5 ) then
605 
606  end if
607  else if ( p%EFTflag == 2 ) then ! designer EFT models
608  if ( p%DesignerEFTmodel == 1 ) then
609  num_param = num_param + 1 ! B0
610  if ( num_param == param_number ) then; param_name = 'B0' ; if ( present(param_name_latex) ) param_name_latex = 'B_0'; return ;
611  end if
612  else if ( p%DesignerEFTmodel == 2 ) then
613  num_param = num_param + 0 ! mc 5e models
614  end if
615  else if ( p%EFTflag == 3 ) then ! EFT alternative parametrizations
616  if ( p%AltParEFTmodel == 1 ) then
617  if ( p%RPHmassPmodel == 1 ) then
618  num_param = num_param + 1 ! RPHmassP0
619  if ( num_param == param_number ) then; param_name = 'RPHmassP0' ; if ( present(param_name_latex) ) param_name_latex = '\tilde{M_0}'; return ;
620  end if
621  else if ( p%RPHmassPmodel == 2 ) then
622  num_param = num_param + 1 ! RPHmassP0
623  if ( num_param == param_number ) then; param_name = 'RPHmassP0' ; if ( present(param_name_latex) ) param_name_latex = '\tilde{M_0}'; return ;
624  end if
625  num_param = num_param + 1 ! RPHmassPexp
626  if ( num_param == param_number ) then; param_name = 'RPHmassPexp' ; if ( present(param_name_latex) ) param_name_latex = '\tilde{M}_{\rm exp}^{\rm RPH}'; return ;
627  end if
628  end if
629  if ( p%RPHkineticitymodel == 1 ) then
630  num_param = num_param + 1 ! RPHkineticity0
631  if ( num_param == param_number ) then; param_name = 'RPHkineticity0' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_0^{\rm K}'; return ;
632  end if
633  else if ( p%RPHkineticitymodel == 2 ) then
634  num_param = num_param + 1 ! RPHkineticity0
635  if ( num_param == param_number ) then; param_name = 'RPHkineticity0' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_0^{\rm K}'; return ;
636  end if
637  num_param = num_param + 1 ! RPHkineticityexp
638  if ( num_param == param_number ) then; param_name = 'RPHkineticityexp' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_{\rm exp}^{\rm K}'; return ;
639  end if
640  end if
641  if ( p%RPHbraidingmodel == 1 ) then
642  num_param = num_param + 1 ! RPHbraiding0
643  if ( num_param == param_number ) then; param_name = 'RPHbraiding0' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_0^{\rm B}'; return ;
644  end if
645  else if ( p%RPHbraidingmodel == 2 ) then
646  num_param = num_param + 1 ! RPHbraiding0
647  if ( num_param == param_number ) then; param_name = 'RPHbraiding0' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_0^{\rm B}'; return ;
648  end if
649  num_param = num_param + 1 ! RPHbraidingexp
650  if ( num_param == param_number ) then; param_name = 'RPHbraidingexp' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_{\rm exp}^{\rm B}'; return ;
651  end if
652  end if
653  if ( p%RPHtensormodel == 1 ) then
654  num_param = num_param + 1 ! RPHtensor0
655  if ( num_param == param_number ) then; param_name = 'RPHtensor0' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_0^{\rm T}'; return ;
656  end if
657  else if ( p%RPHtensormodel == 2 ) then
658  num_param = num_param + 1 ! RPHtensor0
659  if ( num_param == param_number ) then; param_name = 'RPHtensor0' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_0^{\rm T}'; return ;
660  end if
661  num_param = num_param + 1 ! RPHtensorexp
662  if ( num_param == param_number ) then; param_name = 'RPHtensorexp' ; if ( present(param_name_latex) ) param_name_latex = '\alpha_{\rm exp}^{\rm T}'; return ;
663  end if
664  end if
665  end if
666  else if ( p%EFTflag == 4 ) then ! Full mapping EFT model
667  if ( p%FullMappingEFTmodel == 1 ) then
668  if ( p%HoravaSolarSystem .eqv. .true. ) then
669  num_param = num_param + 1 ! Horava_lambda
670  if ( num_param == param_number ) then; param_name = 'Horava_lambda' ; if ( present(param_name_latex) ) param_name_latex = '{\rm Horava}_{\lambda}'; return ;
671  end if
672  num_param = num_param + 1 ! Horava_eta
673  if ( num_param == param_number ) then; param_name = 'Horava_eta' ; if ( present(param_name_latex) ) param_name_latex = '{\rm Horava}_{\eta}'; return ;
674  end if
675  else
676  num_param = num_param + 1 ! Horava_lambda
677  if ( num_param == param_number ) then; param_name = 'Horava_lambda' ; if ( present(param_name_latex) ) param_name_latex = '{\rm Horava}_{\lambda}'; return ;
678  end if
679  num_param = num_param + 1 ! Horava_eta
680  if ( num_param == param_number ) then; param_name = 'Horava_eta' ; if ( present(param_name_latex) ) param_name_latex = '{\rm Horava}_{\eta}'; return ;
681  end if
682  num_param = num_param + 1 ! Horava_xi
683  if ( num_param == param_number ) then; param_name = 'Horava_xi' ; if ( present(param_name_latex) ) param_name_latex = '{\rm Horava}_{\xi}'; return ;
684  end if
685  end if
686  end if
687  end if
688 #endif
689 
690 #ifdef COSMICFISH_MGCAMB
691  if ( p%MGC_model == 1 ) then
692  num_param = num_param + 1 ! B1
693  if ( num_param == param_number ) then; param_name = 'B1' ; if ( present(param_name_latex) ) param_name_latex = 'B_{1}'; return ;
694  end if
695  num_param = num_param + 1 ! B2
696  if ( num_param == param_number ) then; param_name = 'B2' ; if ( present(param_name_latex) ) param_name_latex = 'B_{2}'; return ;
697  end if
698  num_param = num_param + 1 ! lambda1_2
699  if ( num_param == param_number ) then; param_name = 'lambda1_2' ; if ( present(param_name_latex) ) param_name_latex = '\lambda_{1}^{2}'; return ;
700  end if
701  num_param = num_param + 1 ! lambda2_2
702  if ( num_param == param_number ) then; param_name = 'lambda2_2' ; if ( present(param_name_latex) ) param_name_latex = '\lambda_{2}^{2}'; return ;
703  end if
704  num_param = num_param + 1 ! ss
705  if ( num_param == param_number ) then; param_name = 'ss' ; if ( present(param_name_latex) ) param_name_latex = 's'; return ;
706  end if
707  else if ( p%MGC_model == 2 ) then
708  num_param = num_param + 1 ! MGQfix
709  if ( num_param == param_number ) then; param_name = 'MGQfix' ; if ( present(param_name_latex) ) param_name_latex = 'Q_{0}'; return ;
710  end if
711  num_param = num_param + 1 ! MGRfix
712  if ( num_param == param_number ) then; param_name = 'MGRfix' ; if ( present(param_name_latex) ) param_name_latex = 'R_{0}'; return ;
713  end if
714  else if ( p%MGC_model == 3 ) then
715  num_param = num_param + 1 ! Qnot
716  if ( num_param == param_number ) then; param_name = 'Qnot' ; if ( present(param_name_latex) ) param_name_latex = 'Q_{0}'; return ;
717  end if
718  num_param = num_param + 1 ! Rnot
719  if ( num_param == param_number ) then; param_name = 'Rnot' ; if ( present(param_name_latex) ) param_name_latex = 'R_{0}'; return ;
720  end if
721  num_param = num_param + 1 ! sss
722  if ( num_param == param_number ) then; param_name = 'sss' ; if ( present(param_name_latex) ) param_name_latex = 's'; return ;
723  end if
724  else if ( p%MGC_model == 4 ) then
725  num_param = num_param + 1 ! B0
726  if ( num_param == param_number ) then; param_name = 'B0' ; if ( present(param_name_latex) ) param_name_latex = 'B_{0}'; return ;
727  end if
728  else if ( p%MGC_model ==5 ) then
729  num_param = num_param + 1 ! B0
730  if ( num_param == param_number ) then; param_name = 'B0' ; if ( present(param_name_latex) ) param_name_latex = 'B_{0}'; return ;
731  end if
732  num_param = num_param + 1 ! beta1
733  if ( num_param == param_number ) then; param_name = 'beta1' ; if ( present(param_name_latex) ) param_name_latex = '\beta_{1}'; return ;
734  end if
735  num_param = num_param + 1 ! s
736  if ( num_param == param_number ) then; param_name = 's' ; if ( present(param_name_latex) ) param_name_latex = 's'; return ;
737  end if
738  else if ( p%MGC_model ==6 ) then
739  num_param = num_param + 1 ! Linder_gamma
740  if ( num_param == param_number ) then; param_name = 'Linder_gamma' ; if ( present(param_name_latex) ) param_name_latex = '\gamma_{\rm L}'; return ;
741  end if
742  else if ( p%MGC_model == 7 ) then
743  num_param = num_param + 1 ! beta_star
744  if ( num_param == param_number ) then; param_name = 'beta_star' ; if ( present(param_name_latex) ) param_name_latex = '\beta_{\star}'; return ;
745  end if
746  num_param = num_param + 1 ! xi_star
747  if ( num_param == param_number ) then; param_name = 'xi_star' ; if ( present(param_name_latex) ) param_name_latex = '\xi_{\star}'; return ;
748  end if
749  num_param = num_param + 1 ! a_star
750  if ( num_param == param_number ) then; param_name = 'a_star' ; if ( present(param_name_latex) ) param_name_latex = 'a_{\star}'; return ;
751  end if
752  else if ( p%MGC_model == 8 ) then
753  num_param = num_param + 1 ! beta0
754  if ( num_param == param_number ) then; param_name = 'beta0' ; if ( present(param_name_latex) ) param_name_latex = '\beta_{0}'; return ;
755  end if
756  num_param = num_param + 1 ! xi0
757  if ( num_param == param_number ) then; param_name = 'xi0' ; if ( present(param_name_latex) ) param_name_latex = '\xi_{0}'; return ;
758  end if
759  num_param = num_param + 1 ! DilR
760  if ( num_param == param_number ) then; param_name = 'DilR' ; if ( present(param_name_latex) ) param_name_latex = 'R_{\rm dil}'; return ;
761  end if
762  num_param = num_param + 1 ! DilS
763  if ( num_param == param_number ) then; param_name = 'DilS' ; if ( present(param_name_latex) ) param_name_latex = 'S_{\rm dil}'; return ;
764  end if
765  else if ( p%MGC_model == 9 ) then
766  num_param = num_param + 1 ! F_R0
767  if ( num_param == param_number ) then; param_name = 'F_R0' ; if ( present(param_name_latex) ) param_name_latex = 'f_{R,0}'; return ;
768  end if
769  num_param = num_param + 1 ! FRn
770  if ( num_param == param_number ) then; param_name = 'FRn' ; if ( present(param_name_latex) ) param_name_latex = 'n'; return ;
771  end if
772  else if ( p%MGC_model ==10 ) then
773  num_param = num_param + 1 ! beta0
774  if ( num_param == param_number ) then; param_name = 'beta0' ; if ( present(param_name_latex) ) param_name_latex = '\beta_0'; return ;
775  end if
776  num_param = num_param + 1 ! A2
777  if ( num_param == param_number ) then; param_name = 'A2' ; if ( present(param_name_latex) ) param_name_latex = 'A_2'; return ;
778  end if
779  end if
780 #endif
781 
782  end subroutine fisher_param_names
783 
784  ! ---------------------------------------------------------------------------------------------
786  subroutine camb_params_to_params_array( P, FP, num_param, params_array )
788  implicit none
789 
790  Type(cambparams) , intent(in) :: P
791  Type(cosmicfish_params) , intent(in) :: FP
792  integer , intent(in) :: num_param
793  real(dl), dimension(num_param) , intent(out) :: params_array
794 
795  integer :: ind, ind2, i, j, k
796 
797  ind = 1
798 
799  ! standard parameters:
800  ! omegab
801  if ( fp%fisher_par%want_ombh2 ) then
802  params_array(ind) = p%omegab*(p%h0/100._dl)**2
803  ind = ind + 1
804  end if
805  ! omegac
806  if ( fp%fisher_par%want_omch2 ) then
807  params_array(ind) = p%omegac*(p%h0/100._dl)**2
808  ind = ind + 1
809  end if
810  ! omegan
811  if ( fp%fisher_par%want_omnuh2 ) then
812  params_array(ind) = p%omegan*(p%h0/100._dl)**2
813  ind = ind + 1
814  end if
815  ! hubble
816  if ( fp%fisher_par%want_hubble ) then
817  params_array(ind) = p%h0/100._dl
818  ind = ind + 1
819  end if
820  ! helium abundance
821  if ( fp%fisher_par%want_helium_fraction ) then
822  params_array(ind) = p%yhe
823  ind = ind + 1
824  end if
825  ! massless neutrinos
826  if ( fp%fisher_par%want_massless ) then
827  params_array(ind) = p%Num_Nu_massless
828  ind = ind + 1
829  end if
830  ! primordial spectrum parameters:
831  ! scalar_amp
832  if ( fp%fisher_par%want_scalar_amp ) then
833  params_array(ind) = log(10._dl**10) + log(p%InitPower%ScalarPowerAmp(1))
834  ind = ind + 1
835  end if
836  ! scalar_spectral_index
837  if ( fp%fisher_par%want_scalar_spectral_index ) then
838  params_array(ind) = p%InitPower%an(1)
839  ind = ind + 1
840  end if
841  ! scalar_nrun
842  if ( fp%fisher_par%want_scalar_nrun ) then
843  params_array(ind) = p%InitPower%n_run(1)
844  ind = ind + 1
845  end if
846  ! tensor_spectral_index
847  if ( fp%fisher_par%want_tensor_spectral_index ) then
848  params_array(ind) = p%InitPower%ant(1)
849  ind = ind + 1
850  end if
851  ! initial_ratio
852  if ( fp%fisher_par%want_initial_ratio ) then
853  params_array(ind) = p%InitPower%rat(1)
854  ind = ind + 1
855  end if
856  ! reionization: reionization optical depth
857  if ( fp%fisher_par%want_re_optical_depth ) then
858  params_array(ind) = p%Reion%optical_depth
859  ind = ind + 1
860  end if
861 
862  ! Nuisance parameters for LSS:
863  ! bias in every counts redshift windows:
864  if ( fp%fisher_par%want_bias ) then
865  ind2 = 0
866  do i = 1, num_redshiftwindows
867  if ( redshift_w(i)%kind == window_counts ) then
868  ind2 = ind2 +1
869  params_array(ind) = fp%fisher_cls%bias(ind2)
870  ind = ind +1
871  end if
872  end do
873  end if
874 
875  ! nuissance parameters for supernovae:
876  ! alpha_SN
877  if ( fp%fisher_par%want_alpha_SN ) then
878  params_array(ind) = fp%fisher_SN%alpha_SN
879  ind = ind +1
880  end if
881  ! beta_SN
882  if ( fp%fisher_par%want_beta_SN ) then
883  params_array(ind) = fp%fisher_SN%beta_SN
884  ind = ind +1
885  end if
886  ! M0_SN
887  if ( fp%fisher_par%want_M0_SN ) then
888  params_array(ind) = fp%fisher_SN%M0_SN
889  ind = ind +1
890  end if
891 
892 #ifdef COSMICFISH_CAMB
893  ! ppf parameters:
894  ! w0_ppf
895  if ( fp%fisher_par%want_w0_ppf ) then
896  params_array(ind) = p%w_lam
897  ind = ind +1
898  end if
899  ! wa_ppf
900  if ( fp%fisher_par%want_wa_ppf ) then
901  params_array(ind) = p%wa_ppf
902  ind = ind +1
903  end if
904  ! cs_ppf
905  if ( fp%fisher_par%want_cs_ppf ) then
906  params_array(ind) = p%cs2_lam
907  ind = ind +1
908  end if
909 #endif
910 
911 #ifdef COSMICFISH_EFTCAMB
912  ! EFTCAMB parameters
913  ! w_DE parameters:
914  if ( p%EFTflag /= 0 .and. p%EFTwDE /= 0 ) then
915  if ( p%EFTwDE == 1 ) then ! DE with constant Eos
916  params_array(ind) = p%EFTw0
917  ind = ind +1
918  else if ( p%EFTwDE == 2 ) then ! CPL parametrization
919  params_array(ind) = p%EFTw0
920  params_array(ind+1) = p%EFTwa
921  ind = ind +2
922  else if ( p%EFTwDE == 3 ) then ! JBP parametrization
923  params_array(ind) = p%EFTw0
924  params_array(ind+1) = p%EFTwa
925  params_array(ind+2) = p%EFTwn
926  ind = ind +3
927  else if ( p%EFTwDE == 4 ) then ! turning point parametrization
928  params_array(ind) = p%EFTw0
929  params_array(ind+1) = p%EFTwa
930  params_array(ind+2) = p%EFTwat
931  ind = ind +3
932  else if ( p%EFTwDE == 5 ) then ! Taylor expansion
933  params_array(ind) = p%EFTw0
934  params_array(ind+1) = p%EFTwa
935  params_array(ind+2) = p%EFTw2
936  params_array(ind+3) = p%EFTw3
937  ind = ind +4
938  else if ( p%EFTwDE == 6 ) then ! User defined
939  stop 'You have to define the parametrization first'
940  end if
941  end if
942 
943  ! EFT functions parameters:
944  if ( p%EFTflag == 1 ) then ! Pure EFT models
945  ! EFTOmega:
946  if ( p%PureEFTmodelOmega == 1 ) then
947  params_array(ind) = p%EFTOmega0
948  ind = ind +1
949  else if ( p%PureEFTmodelOmega == 2 ) then
950  params_array(ind) = p%EFTOmega0
951  ind = ind +1
952  else if ( p%PureEFTmodelOmega == 3 ) then
953  params_array(ind) = p%EFTOmega0
954  params_array(ind+1) = p%EFTOmegaExp
955  ind = ind +2
956  else if ( p%PureEFTmodelOmega == 4 ) then
957  params_array(ind) = p%EFTOmega0
958  params_array(ind+1) = p%EFTOmegaExp
959  ind = ind +2
960  else if ( p%PureEFTmodelOmega == 5 ) then
961  stop 'You have to define the parametrization first'
962  end if
963  ! EFTGamma1
964  if ( p%PureEFTmodelGamma1 == 1 ) then
965  params_array(ind) = p%EFTGamma10
966  ind = ind +1
967  else if ( p%PureEFTmodelGamma1 == 2 ) then
968  params_array(ind) = p%EFTGamma10
969  ind = ind +1
970  else if ( p%PureEFTmodelGamma1 == 3 ) then
971  params_array(ind) = p%EFTGamma10
972  params_array(ind+1) = p%EFTGamma1Exp
973  ind = ind +2
974  else if ( p%PureEFTmodelGamma1 == 4 ) then
975  params_array(ind) = p%EFTGamma10
976  params_array(ind+1) = p%EFTGamma1Exp
977  ind = ind +2
978  else if ( p%PureEFTmodelGamma1 == 5 ) then
979  stop 'You have to define the parametrization first'
980  end if
981  ! EFTGamma2
982  if ( p%PureEFTmodelGamma2 == 1 ) then
983  params_array(ind) = p%EFTGamma20
984  ind = ind +1
985  else if ( p%PureEFTmodelGamma2 == 2 ) then
986  params_array(ind) = p%EFTGamma20
987  ind = ind +1
988  else if ( p%PureEFTmodelGamma2 == 3 ) then
989  params_array(ind) = p%EFTGamma20
990  params_array(ind+1) = p%EFTGamma2Exp
991  ind = ind +2
992  else if ( p%PureEFTmodelGamma2 == 4 ) then
993  params_array(ind) = p%EFTGamma20
994  params_array(ind+1) = p%EFTGamma2Exp
995  ind = ind +2
996  else if ( p%PureEFTmodelGamma2 == 5 ) then
997  stop 'You have to define the parametrization first'
998  end if
999  ! EFTGamma3
1000  if ( p%PureEFTmodelGamma3 == 1 ) then
1001  params_array(ind) = p%EFTGamma30
1002  ind = ind +1
1003  else if ( p%PureEFTmodelGamma3 == 2 ) then
1004  params_array(ind) = p%EFTGamma30
1005  ind = ind +1
1006  else if ( p%PureEFTmodelGamma3 == 3 ) then
1007  params_array(ind) = p%EFTGamma30
1008  params_array(ind+1) = p%EFTGamma3Exp
1009  ind = ind +2
1010  else if ( p%PureEFTmodelGamma3 == 4 ) then
1011  params_array(ind) = p%EFTGamma30
1012  params_array(ind+1) = p%EFTGamma3Exp
1013  ind = ind +2
1014  else if ( p%PureEFTmodelGamma3 == 5 ) then
1015  stop 'You have to define the parametrization first'
1016  end if
1017  ! EFTGamma4
1018  if ( p%PureEFTmodelGamma4 == 1 ) then
1019  params_array(ind) = p%EFTGamma40
1020  ind = ind +1
1021  else if ( p%PureEFTmodelGamma4 == 2 ) then
1022  params_array(ind) = p%EFTGamma40
1023  ind = ind +1
1024  else if ( p%PureEFTmodelGamma4 == 3 ) then
1025  params_array(ind) = p%EFTGamma40
1026  params_array(ind+1) = p%EFTGamma4Exp
1027  ind = ind +2
1028  else if ( p%PureEFTmodelGamma4 == 4 ) then
1029  params_array(ind) = p%EFTGamma40
1030  params_array(ind+1) = p%EFTGamma4Exp
1031  ind = ind +2
1032  else if ( p%PureEFTmodelGamma4 == 5 ) then
1033  stop 'You have to define the parametrization first'
1034  end if
1035  ! EFTGamma5
1036  if ( p%PureEFTmodelGamma5 == 1 ) then
1037  params_array(ind) = p%EFTGamma50
1038  ind = ind +1
1039  else if ( p%PureEFTmodelGamma5 == 2 ) then
1040  params_array(ind) = p%EFTGamma50
1041  ind = ind +1
1042  else if ( p%PureEFTmodelGamma5 == 3 ) then
1043  params_array(ind) = p%EFTGamma50
1044  params_array(ind+1) = p%EFTGamma5Exp
1045  ind = ind +2
1046  else if ( p%PureEFTmodelGamma5 == 4 ) then
1047  params_array(ind) = p%EFTGamma50
1048  params_array(ind+1) = p%EFTGamma5Exp
1049  ind = ind +2
1050  else if ( p%PureEFTmodelGamma5 == 5 ) then
1051  stop 'You have to define the parametrization first'
1052  end if
1053  ! EFTGamma6
1054  if ( p%PureEFTmodelGamma6 == 1 ) then
1055  params_array(ind) = p%EFTGamma60
1056  ind = ind +1
1057  else if ( p%PureEFTmodelGamma6 == 2 ) then
1058  params_array(ind) = p%EFTGamma60
1059  ind = ind +1
1060  else if ( p%PureEFTmodelGamma6 == 3 ) then
1061  params_array(ind) = p%EFTGamma60
1062  params_array(ind+1) = p%EFTGamma6Exp
1063  ind = ind +2
1064  else if ( p%PureEFTmodelGamma6 == 4 ) then
1065  params_array(ind) = p%EFTGamma60
1066  params_array(ind+1) = p%EFTGamma6Exp
1067  ind = ind +2
1068  else if ( p%PureEFTmodelGamma6 == 5 ) then
1069  stop 'You have to define the parametrization first'
1070  end if
1071  else if ( p%EFTflag == 2 ) then ! designer EFT models
1072  if ( p%DesignerEFTmodel == 1 ) then
1073  params_array(ind) = p%EFTB0
1074  ind = ind +1
1075  else if ( p%DesignerEFTmodel == 2 ) then
1076  end if
1077  else if ( p%EFTflag == 3 ) then ! EFT alternative parametrizations
1078  if ( p%AltParEFTmodel == 1 ) then
1079  if ( p%RPHmassPmodel == 1 ) then
1080  params_array(ind) = p%RPHmassP0
1081  ind = ind +1
1082  else if ( p%RPHmassPmodel == 2 ) then
1083  params_array(ind) = p%RPHmassP0
1084  params_array(ind+1) = p%RPHmassPexp
1085  ind = ind +2
1086  end if
1087  if ( p%RPHkineticitymodel == 1 ) then
1088  params_array(ind) = p%RPHkineticity0
1089  ind = ind +1
1090  else if ( p%RPHkineticitymodel == 2 ) then
1091  params_array(ind) = p%RPHkineticity0
1092  params_array(ind+1) = p%RPHkineticityexp
1093  ind = ind +2
1094  end if
1095  if ( p%RPHbraidingmodel == 1 ) then
1096  params_array(ind) = p%RPHbraiding0
1097  ind = ind +1
1098  else if ( p%RPHbraidingmodel == 2 ) then
1099  params_array(ind) = p%RPHbraiding0
1100  params_array(ind+1) = p%RPHbraidingexp
1101  ind = ind +2
1102  end if
1103  if ( p%RPHtensormodel == 1 ) then
1104  params_array(ind) = p%RPHtensor0
1105  ind = ind +1
1106  else if ( p%RPHtensormodel == 2 ) then
1107  params_array(ind) = p%RPHtensor0
1108  params_array(ind+1) = p%RPHtensorexp
1109  ind = ind +2
1110  end if
1111  end if
1112  else if ( p%EFTflag == 4 ) then ! Full mapping EFT model
1113  if ( p%FullMappingEFTmodel == 1 ) then
1114  if ( p%HoravaSolarSystem .eqv. .true. ) then
1115  params_array(ind) = p%Horava_lambda
1116  params_array(ind+1) = p%Horava_eta
1117  ind = ind +2
1118  else
1119  params_array(ind) = p%Horava_lambda
1120  params_array(ind+1) = p%Horava_eta
1121  params_array(ind+2) = p%Horava_xi
1122  ind = ind +3
1123  end if
1124  end if
1125  end if
1126 #endif
1127 
1128 #ifdef COSMICFISH_MGCAMB
1129  if ( p%MGC_model == 1 ) then
1130  params_array(ind) = p%B1
1131  params_array(ind+1) = p%B2
1132  params_array(ind+2) = p%lambda1_2
1133  params_array(ind+3) = p%lambda2_2
1134  params_array(ind+4) = p%ss
1135  ind = ind +5
1136  else if ( p%MGC_model == 2 ) then
1137  params_array(ind) = p%MGQfix
1138  params_array(ind+1) = p%MGRfix
1139  ind = ind +2
1140  else if ( p%MGC_model == 3 ) then
1141  params_array(ind) = p%Qnot
1142  params_array(ind+1) = p%Rnot
1143  params_array(ind+2) = p%sss
1144  ind = ind +3
1145  else if ( p%MGC_model == 4 ) then
1146  params_array(ind) = p%B0
1147  ind = ind +1
1148  else if ( p%MGC_model ==5 ) then
1149  params_array(ind) = p%B0
1150  params_array(ind+1) = p%B1
1151  params_array(ind+2) = p%ss
1152  ind = ind +3
1153  else if ( p%MGC_model ==6 ) then
1154  params_array(ind) = p%Linder_gamma
1155  ind = ind +1
1156  else if ( p%MGC_model == 7 ) then
1157  params_array(ind) = p%beta_star
1158  params_array(ind+1) = p%xi_star
1159  params_array(ind+2) = p%a_star
1160  ind = ind +3
1161  else if ( p%MGC_model == 8 ) then
1162  params_array(ind) = p%beta0
1163  params_array(ind+1) = p%xi0
1164  params_array(ind+2) = p%DilR
1165  params_array(ind+3) = p%DilS
1166  ind = ind +4
1167  else if ( p%MGC_model == 9 ) then
1168  params_array(ind) = p%F_R0
1169  params_array(ind+1) = p%FRn
1170  ind = ind +2
1171  else if ( p%MGC_model ==10 ) then
1172  params_array(ind) = p%beta0
1173  params_array(ind+1) = p%A_2
1174  ind = ind +2
1175  end if
1176 #endif
1177 
1178  ! Quality check:
1179  if ( ind-1 /= num_param ) stop 'CAMB_params_to_params_array: Wrong assignment of parameters'
1180 
1181  end subroutine camb_params_to_params_array
1182 
1183  ! ---------------------------------------------------------------------------------------------
1186  subroutine params_array_to_camb_param( P, FP, num_param, params_array )
1188  implicit none
1189 
1190  Type(cambparams) :: P
1191  Type(cosmicfish_params) :: FP
1192  integer , intent(in) :: num_param
1193  real(dl), dimension(num_param) , intent(in) :: params_array
1194 
1195  integer :: ind, ind2, i, j, k
1196  real(dl) :: temp
1197 
1198  ind = 1
1199  i = 1
1200 
1201  if ( fp%fisher_par%want_hubble ) then
1202  if ( fp%fisher_par%want_ombh2 ) i = i + 1
1203  if ( fp%fisher_par%want_omch2 ) i = i + 1
1204  if ( fp%fisher_par%want_omnuh2 ) i = i + 1
1205  p%h0 = params_array(i)*100._dl
1206  end if
1207 
1208  if ( fp%fisher_par%want_ombh2 ) then
1209  p%omegab = params_array(ind)/(p%h0/100._dl)**2
1210  ind = ind + 1
1211  end if
1212 
1213  if ( fp%fisher_par%want_omch2 ) then
1214  p%omegac = params_array(ind)/(p%h0/100._dl)**2
1215  ind = ind + 1
1216  end if
1217 
1218  if ( fp%fisher_par%want_omnuh2 ) then
1219  p%omegan = params_array(ind)/(p%h0/100._dl)**2
1220  ind = ind + 1
1221  end if
1222 
1223  if ( fp%fisher_par%want_hubble ) then
1224  ind = ind + 1
1225  end if
1226 
1227  p%omegav = 1- p%omegab- p%omegac - p%omegan
1228 
1229  ! helium abundance
1230  if ( fp%fisher_par%want_helium_fraction ) then
1231  p%yhe = params_array(ind)
1232  ind = ind + 1
1233  end if
1234  ! massless neutrinos
1235  if ( fp%fisher_par%want_massless ) then
1236  p%Num_Nu_massless = params_array(ind)
1237  ind = ind + 1
1238  end if
1239  ! primordial spectrum parameters:
1240  ! scalar_amp
1241  if ( fp%fisher_par%want_scalar_amp ) then
1242  p%InitPower%ScalarPowerAmp(1) = exp(params_array(ind))/(10._dl**10)
1243  ind = ind + 1
1244  end if
1245  ! scalar_spectral_index
1246  if ( fp%fisher_par%want_scalar_spectral_index ) then
1247  p%InitPower%an(1) = params_array(ind)
1248  ind = ind + 1
1249  end if
1250  ! scalar_nrun
1251  if ( fp%fisher_par%want_scalar_nrun ) then
1252  p%InitPower%n_run(1) = params_array(ind)
1253  ind = ind + 1
1254  end if
1255  ! tensor_spectral_index
1256  if ( fp%fisher_par%want_tensor_spectral_index ) then
1257  p%InitPower%ant(1) = params_array(ind)
1258  ind = ind + 1
1259  end if
1260  ! initial_ratio
1261  if ( fp%fisher_par%want_initial_ratio ) then
1262  p%InitPower%rat(1) = params_array(ind)
1263  ind = ind + 1
1264  end if
1265  ! reionization: reionization optical depth
1266  if ( fp%fisher_par%want_re_optical_depth ) then
1267  p%Reion%optical_depth = params_array(ind)
1268  ind = ind + 1
1269  end if
1270 
1271  ! Nuisance parameters for LSS:
1272  ! bias in every counts redshift windows:
1273  if ( fp%fisher_par%want_bias ) then
1274  ind2 = 0
1275  do i = 1, num_redshiftwindows
1276  if ( redshift_w(i)%kind == window_counts ) then
1277  ind2 = ind2 +1
1278  fp%fisher_cls%bias(ind2) = params_array(ind)
1279  redshift_w(i)%bias = fp%fisher_cls%bias(ind2)
1280  ind = ind +1
1281  end if
1282  end do
1283  end if
1284 
1285  ! nuissance parameters for supernovae:
1286  ! alpha_SN
1287  if ( fp%fisher_par%want_alpha_SN ) then
1288  fp%fisher_SN%alpha_SN = params_array(ind)
1289  ind = ind +1
1290  end if
1291  ! beta_SN
1292  if ( fp%fisher_par%want_beta_SN ) then
1293  fp%fisher_SN%beta_SN = params_array(ind)
1294  ind = ind +1
1295  end if
1296  ! M0_SN
1297  if ( fp%fisher_par%want_M0_SN ) then
1298  fp%fisher_SN%M0_SN = params_array(ind)
1299  ind = ind +1
1300  end if
1301 
1302 #ifdef COSMICFISH_CAMB
1303  ! ppf parameters:
1304  ! w0_ppf
1305  if ( fp%fisher_par%want_w0_ppf ) then
1306  p%w_lam = params_array(ind)
1307  ind = ind +1
1308  end if
1309  ! wa_ppf
1310  if ( fp%fisher_par%want_wa_ppf ) then
1311  p%wa_ppf = params_array(ind)
1312  ind = ind +1
1313  end if
1314  ! cs_ppf
1315  if ( fp%fisher_par%want_cs_ppf ) then
1316  p%cs2_lam = params_array(ind)
1317  call setcgammappf( p )
1318  ind = ind +1
1319  end if
1320 #endif
1321 
1322 #ifdef COSMICFISH_EFTCAMB
1323  ! EFTCAMB parameters
1324  ! w_DE parameters:
1325  if ( p%EFTflag /= 0 .and. p%EFTwDE /= 0 ) then
1326  if ( p%EFTwDE == 1 ) then ! DE with constant Eos
1327  p%EFTw0 = params_array(ind)
1328  ind = ind +1
1329  else if ( p%EFTwDE == 2 ) then ! CPL parametrization
1330  p%EFTw0 = params_array(ind)
1331  p%EFTwa = params_array(ind+1)
1332  ind = ind +2
1333  else if ( p%EFTwDE == 3 ) then ! JBP parametrization
1334  p%EFTw0 = params_array(ind)
1335  p%EFTwa = params_array(ind+1)
1336  p%EFTwn = params_array(ind+2)
1337  ind = ind +3
1338  else if ( p%EFTwDE == 4 ) then ! turning point parametrization
1339  p%EFTw0 = params_array(ind)
1340  p%EFTwa = params_array(ind+1)
1341  p%EFTwat = params_array(ind+2)
1342  ind = ind +3
1343  else if ( p%EFTwDE == 5 ) then ! Taylor expansion
1344  p%EFTw0 = params_array(ind)
1345  p%EFTwa = params_array(ind+1)
1346  p%EFTw2 = params_array(ind+2)
1347  p%EFTw3 = params_array(ind+3)
1348  ind = ind +4
1349  else if ( p%EFTwDE == 6 ) then ! User defined
1350  stop 'You have to define the parametrization first'
1351  end if
1352  end if
1353 
1354  ! EFT functions parameters:
1355  if ( p%EFTflag == 1 ) then ! Pure EFT models
1356  ! EFTOmega:
1357  if ( p%PureEFTmodelOmega == 1 ) then
1358  p%EFTOmega0 = params_array(ind)
1359  ind = ind +1
1360  else if ( p%PureEFTmodelOmega == 2 ) then
1361  p%EFTOmega0 = params_array(ind)
1362  ind = ind +1
1363  else if ( p%PureEFTmodelOmega == 3 ) then
1364  p%EFTOmega0 = params_array(ind)
1365  p%EFTOmegaExp = params_array(ind+1)
1366  ind = ind +2
1367  else if ( p%PureEFTmodelOmega == 4 ) then
1368  p%EFTOmega0 = params_array(ind)
1369  p%EFTOmegaExp = params_array(ind+1)
1370  ind = ind +2
1371  else if ( p%PureEFTmodelOmega == 5 ) then
1372  stop 'You have to define the parametrization first'
1373  end if
1374  ! EFTGamma1
1375  if ( p%PureEFTmodelGamma1 == 1 ) then
1376  p%EFTGamma10 = params_array(ind)
1377  ind = ind +1
1378  else if ( p%PureEFTmodelGamma1 == 2 ) then
1379  p%EFTGamma10 = params_array(ind)
1380  ind = ind +1
1381  else if ( p%PureEFTmodelGamma1 == 3 ) then
1382  p%EFTGamma10 = params_array(ind)
1383  p%EFTGamma1Exp = params_array(ind+1)
1384  ind = ind +2
1385  else if ( p%PureEFTmodelGamma1 == 4 ) then
1386  p%EFTGamma10 = params_array(ind)
1387  p%EFTGamma1Exp = params_array(ind+1)
1388  ind = ind +2
1389  else if ( p%PureEFTmodelGamma1 == 5 ) then
1390  stop 'You have to define the parametrization first'
1391  end if
1392  ! EFTGamma2
1393  if ( p%PureEFTmodelGamma2 == 1 ) then
1394  p%EFTGamma20 = params_array(ind)
1395  ind = ind +1
1396  else if ( p%PureEFTmodelGamma2 == 2 ) then
1397  p%EFTGamma20 = params_array(ind)
1398  ind = ind +1
1399  else if ( p%PureEFTmodelGamma2 == 3 ) then
1400  p%EFTGamma20 = params_array(ind)
1401  p%EFTGamma2Exp = params_array(ind+1)
1402  ind = ind +2
1403  else if ( p%PureEFTmodelGamma2 == 4 ) then
1404  p%EFTGamma20 = params_array(ind)
1405  p%EFTGamma2Exp = params_array(ind+1)
1406  ind = ind +2
1407  else if ( p%PureEFTmodelGamma2 == 5 ) then
1408  stop 'You have to define the parametrization first'
1409  end if
1410  ! EFTGamma3
1411  if ( p%PureEFTmodelGamma3 == 1 ) then
1412  p%EFTGamma30 = params_array(ind)
1413  ind = ind +1
1414  else if ( p%PureEFTmodelGamma3 == 2 ) then
1415  p%EFTGamma30 = params_array(ind)
1416  ind = ind +1
1417  else if ( p%PureEFTmodelGamma3 == 3 ) then
1418  p%EFTGamma30 = params_array(ind)
1419  p%EFTGamma3Exp = params_array(ind+1)
1420  ind = ind +2
1421  else if ( p%PureEFTmodelGamma3 == 4 ) then
1422  p%EFTGamma30 = params_array(ind)
1423  p%EFTGamma3Exp = params_array(ind+1)
1424  ind = ind +2
1425  else if ( p%PureEFTmodelGamma3 == 5 ) then
1426  stop 'You have to define the parametrization first'
1427  end if
1428  ! EFTGamma4
1429  if ( p%PureEFTmodelGamma4 == 1 ) then
1430  p%EFTGamma40 = params_array(ind)
1431  ind = ind +1
1432  else if ( p%PureEFTmodelGamma4 == 2 ) then
1433  p%EFTGamma40 = params_array(ind)
1434  ind = ind +1
1435  else if ( p%PureEFTmodelGamma4 == 3 ) then
1436  p%EFTGamma40 = params_array(ind)
1437  p%EFTGamma4Exp = params_array(ind+1)
1438  ind = ind +2
1439  else if ( p%PureEFTmodelGamma4 == 4 ) then
1440  p%EFTGamma40 = params_array(ind)
1441  p%EFTGamma4Exp = params_array(ind+1)
1442  ind = ind +2
1443  else if ( p%PureEFTmodelGamma4 == 5 ) then
1444  stop 'You have to define the parametrization first'
1445  end if
1446  ! EFTGamma5
1447  if ( p%PureEFTmodelGamma5 == 1 ) then
1448  p%EFTGamma50 = params_array(ind)
1449  ind = ind +1
1450  else if ( p%PureEFTmodelGamma5 == 2 ) then
1451  p%EFTGamma50 = params_array(ind)
1452  ind = ind +1
1453  else if ( p%PureEFTmodelGamma5 == 3 ) then
1454  p%EFTGamma50 = params_array(ind)
1455  p%EFTGamma5Exp = params_array(ind+1)
1456  ind = ind +2
1457  else if ( p%PureEFTmodelGamma5 == 4 ) then
1458  p%EFTGamma50 = params_array(ind)
1459  p%EFTGamma5Exp = params_array(ind+1)
1460  ind = ind +2
1461  else if ( p%PureEFTmodelGamma5 == 5 ) then
1462  stop 'You have to define the parametrization first'
1463  end if
1464  ! EFTGamma6
1465  if ( p%PureEFTmodelGamma6 == 1 ) then
1466  p%EFTGamma60 = params_array(ind)
1467  ind = ind +1
1468  else if ( p%PureEFTmodelGamma6 == 2 ) then
1469  p%EFTGamma60 = params_array(ind)
1470  ind = ind +1
1471  else if ( p%PureEFTmodelGamma6 == 3 ) then
1472  p%EFTGamma60 = params_array(ind)
1473  p%EFTGamma6Exp = params_array(ind+1)
1474  ind = ind +2
1475  else if ( p%PureEFTmodelGamma6 == 4 ) then
1476  p%EFTGamma60 = params_array(ind)
1477  p%EFTGamma6Exp = params_array(ind+1)
1478  ind = ind +2
1479  else if ( p%PureEFTmodelGamma6 == 5 ) then
1480  stop 'You have to define the parametrization first'
1481  end if
1482  else if ( p%EFTflag == 2 ) then ! designer EFT models
1483  if ( p%DesignerEFTmodel == 1 ) then
1484  p%EFTB0 = params_array(ind)
1485  ind = ind +1
1486  else if ( p%DesignerEFTmodel == 2 ) then
1487  end if
1488  else if ( p%EFTflag == 3 ) then ! EFT alternative parametrizations
1489  if ( p%AltParEFTmodel == 1 ) then
1490  if ( p%RPHmassPmodel == 1 ) then
1491  p%RPHmassP0 = params_array(ind)
1492  ind = ind +1
1493  else if ( p%RPHmassPmodel == 2 ) then
1494  p%RPHmassP0 = params_array(ind)
1495  p%RPHmassPexp = params_array(ind+1)
1496  ind = ind +2
1497  end if
1498  if ( p%RPHkineticitymodel == 1 ) then
1499  p%RPHkineticity0 = params_array(ind)
1500  ind = ind +1
1501  else if ( p%RPHkineticitymodel == 2 ) then
1502  p%RPHkineticity0 = params_array(ind)
1503  p%RPHkineticityexp = params_array(ind+1)
1504  ind = ind +2
1505  end if
1506  if ( p%RPHbraidingmodel == 1 ) then
1507  p%RPHbraiding0 = params_array(ind)
1508  ind = ind +1
1509  else if ( p%RPHbraidingmodel == 2 ) then
1510  p%RPHbraiding0 = params_array(ind)
1511  p%RPHbraidingexp = params_array(ind+1)
1512  ind = ind +2
1513  end if
1514  if ( p%RPHtensormodel == 1 ) then
1515  p%RPHtensor0 = params_array(ind)
1516  ind = ind +1
1517  else if ( p%RPHtensormodel == 2 ) then
1518  p%RPHtensor0 = params_array(ind)
1519  p%RPHtensorexp = params_array(ind+1)
1520  ind = ind +2
1521  end if
1522  end if
1523  else if ( p%EFTflag == 4 ) then ! Full mapping EFT model
1524  if ( p%FullMappingEFTmodel == 1 ) then
1525  if ( p%HoravaSolarSystem .eqv. .true. ) then
1526  p%Horava_lambda = params_array(ind)
1527  p%Horava_eta = params_array(ind+1)
1528  ind = ind +2
1529  else
1530  p%Horava_lambda = params_array(ind)
1531  p%Horava_eta = params_array(ind+1)
1532  p%Horava_xi = params_array(ind+2)
1533  ind = ind +3
1534  end if
1535  end if
1536  end if
1537 #endif
1538 
1539 #ifdef COSMICFISH_MGCAMB
1540  if ( p%MGC_model == 1 ) then
1541  p%B1 = params_array(ind)
1542  p%B2 = params_array(ind+1)
1543  p%lambda1_2 = params_array(ind+2)
1544  p%lambda2_2 = params_array(ind+3)
1545  p%ss = params_array(ind+4)
1546  ind = ind +5
1547  else if ( p%MGC_model == 2 ) then
1548  p%MGQfix = params_array(ind)
1549  p%MGRfix = params_array(ind+1)
1550  ind = ind +2
1551  else if ( p%MGC_model == 3 ) then
1552  p%Qnot = params_array(ind)
1553  p%Rnot = params_array(ind+1)
1554  p%sss = params_array(ind+2)
1555  ind = ind +3
1556  else if ( p%MGC_model == 4 ) then
1557  p%B0 = params_array(ind)
1558  p%B1 = 4.d0/3.d0
1559  p%lambda1_2 = (p%B0*(299792458.d-3)**2)/(2.d0*p%H0**2)
1560  p%B2 = 0.5d0
1561  p%lambda2_2 = p%B1*p%lambda1_2
1562  p%ss = 4.d0
1563  ind = ind +1
1564  else if ( p%MGC_model ==5 ) then
1565  p%B0 = params_array(ind)
1566  p%B1 = params_array(ind+1)
1567  p%lambda1_2 = (p%B0*(299792458.d-3)**2)/(2.d0*p%H0**2)
1568  p%B2 = 2.d0/p%B1 -1.d0
1569  p%lambda2_2 = p%B1*p%lambda1_2
1570  p%ss = params_array(ind+2)
1571  ind = ind +3
1572  else if ( p%MGC_model ==6 ) then
1573  p%Linder_gamma = params_array(ind)
1574  ind = ind +1
1575  else if ( p%MGC_model == 7 ) then
1576  p%beta_star = params_array(ind)
1577  p%xi_star = params_array(ind+1)
1578  p%a_star = params_array(ind+2)
1579  p%GRtrans = p%a_star
1580  ind = ind +3
1581  else if ( p%MGC_model == 8 ) then
1582  p%beta0 = params_array(ind)
1583  p%xi0 = params_array(ind+1)
1584  p%DilR = params_array(ind+2)
1585  p%DilS = params_array(ind+3)
1586  ind = ind +4
1587  else if ( p%MGC_model == 9 ) then
1588  p%F_R0 = params_array(ind)
1589  p%FRn = params_array(ind+1)
1590  p%beta0 = 1.d0/sqrt(6.d0)
1591  ind = ind +2
1592  else if ( p%MGC_model ==10 ) then
1593  p%beta0 = params_array(ind)
1594  p%A_2 = params_array(ind+1)
1595  ind = ind +2
1596  end if
1597 #endif
1598 
1599  ! Quality check:
1600  if ( ind-1 /= num_param ) stop 'Wrong assignment of parameters'
1601  if (.not. camb_validateparams(p)) stop 'Stopped due to parameter error'
1602 
1603  end subroutine params_array_to_camb_param
1604 
1605  ! ---------------------------------------------------------------------------------------------
1608  subroutine save_fisher_to_file( P, FP, fish_dim, fisher_matrix, filename )
1610  implicit none
1611 
1612  Type(cambparams) :: P
1613  Type(cosmicfish_params) :: FP
1614  integer , intent(in) :: fish_dim
1615  real(dl), intent(in), dimension(fish_dim, fish_dim) :: fisher_matrix
1616  character(len=*), intent(in) :: filename
1617 
1618  integer :: ind
1619  character(LEN=500) :: param_name
1620  character(LEN=500) :: param_name_latex
1621  real(dl), allocatable, dimension(:) :: parameters
1622 
1623  allocate( parameters(fish_dim) )
1624  call camb_params_to_params_array( p, fp, fish_dim, parameters )
1625 
1626  open(unit=666, file=filename, action="write", status="replace")
1627 
1628  ! write the header:
1629  write(666,'(a)') '#'
1630  write(666,'(a)') '# This file contains a Fisher matrix created with the CosmicFish code.'
1631  ! write the parameter names:
1632  write(666,'(a)') '#'
1633  write(666,'(a)') '# The parameters of this Fisher matrix are:'
1634  write(666,'(a)') '#'
1635 
1636  do ind = 1, fish_dim
1637  call fisher_param_names(p, fp, ind, param_name, param_name_latex)
1638  write(666,'(a,I5,a,a,a,a,a,E25.16)') '#', ind, ' ', trim(param_name), ' ', trim(param_name_latex), ' ', parameters(ind)
1639  end do
1640 
1641  write(666,'(a)') '#'
1642  write(666,'(a)')
1643 
1644  ! write the Fisher matrix:
1645  do ind = 1, fish_dim
1646  write(666,'(2x,*(E25.16,2x))') fisher_matrix(ind, :)
1647  end do
1648 
1649  close(666)
1650 
1651  end subroutine save_fisher_to_file
1652 
1653  ! ---------------------------------------------------------------------------------------------
1655  subroutine save_paramnames_to_file( P, FP, fish_dim, filename )
1657  implicit none
1658 
1659  Type(cambparams) :: P
1660  Type(cosmicfish_params) :: FP
1661  integer , intent(in) :: fish_dim
1662  character(len=*), intent(in) :: filename
1663 
1664  integer :: ind
1665  character(LEN=500) :: param_name
1666  character(LEN=500) :: param_name_latex
1667  real(dl), allocatable, dimension(:) :: parameters
1668 
1669  allocate( parameters(fish_dim) )
1670  call camb_params_to_params_array( p, fp, fish_dim, parameters )
1671 
1672  open(unit=666, file=filename, action="write", status="replace")
1673 
1674  ! write the header:
1675 
1676  write(666,'(a)') '#'
1677  write(666,'(a)') '# This file contains the parameter names for a Fisher matrix.'
1678  write(666,'(a)') '#'
1679  ! write the parameter names:
1680 
1681  do ind = 1, fish_dim
1682  call fisher_param_names(p, fp, ind, param_name, param_name_latex)
1683  write(666,'(a,a,a,a,E25.16)') trim(param_name), ' ', trim(param_name_latex), ' ', parameters(ind)
1684  end do
1685 
1686  close(666)
1687 
1688  end subroutine save_paramnames_to_file
1689 
1690 end module cosmicfish_camb_interface
subroutine, public params_array_to_camb_param(P, FP, num_param, params_array)
This subroutine updates the values of the parameters P and FP based on the values contained in the pa...
subroutine, public save_fisher_to_file(P, FP, fish_dim, fisher_matrix, filename)
This subroutine saves the Fisher matrix to file taking care of putting into the file all the relevant...
subroutine, public camb_params_to_params_array(P, FP, num_param, params_array)
This subroutine updates the values of the parameters array based on the values of P and FP...
This module contains the definitions of derived data types used to store the informations about the f...
subroutine, public fisher_param_names(P, FP, param_number, param_name, param_name_latex)
This subroutine returns the name corresponding to a parameter number.
This module contains the interface layer between CosmicFish and CAMB.
subroutine, public num_param_fisher(P, FP, num_param)
This subroutine returns the number of parameters in the Fisher matrix based on the input parameters...
subroutine, public save_paramnames_to_file(P, FP, fish_dim, filename)
This subroutine saves a file containing the parameter names for the Fisher matrix run...
This derived data type contains all the informations that the cosmicfish code needs to run...