source: Sophya/trunk/Cosmo/SimLSS/hu_tf_fit.c@ 3986

Last change on this file since 3986 was 3315, checked in by cmv, 18 years ago

add code fct transfert de W.Hu + prog de comparaison avec pkspectrum.cc cmv 24/08/2007

File size: 14.4 KB
RevLine 
[3315]1/* The following routines implement all of the fitting formulae in
2Eisenstein \& Hu (1997) */
3
4/* There are two sets of routines here. The first set,
5
6 TFfit_hmpc(), TFset_parameters(), and TFfit_onek(),
7
8calculate the transfer function for an arbitrary CDM+baryon universe using
9the fitting formula in Section 3 of the paper. The second set,
10
11 TFsound_horizon_fit(), TFk_peak(), TFnowiggles(), and TFzerobaryon(),
12
13calculate other quantities given in Section 4 of the paper. */
14
15#include <math.h>
16#include <stdio.h>
17#include <stdlib.h>
18
19void TFset_parameters(double omega0hh, double f_baryon, double Tcmb);
20double TFfit_onek(double k, double *tf_baryon, double *tf_cdm);
21
22void TFfit_hmpc(double omega0, double f_baryon, double hubble, double Tcmb,
23 int numk, double *k, double *tf_full, double *tf_baryon, double *tf_cdm);
24
25double TFsound_horizon_fit(double omega0, double f_baryon, double hubble);
26double TFk_peak(double omega0, double f_baryon, double hubble);
27double TFnowiggles(double omega0, double f_baryon, double hubble,
28 double Tcmb, double k_hmpc);
29double TFzerobaryon(double omega0, double hubble, double Tcmb, double k_hmpc);
30
31/* ------------------------ DRIVER ROUTINE --------------------------- */
32/* The following is an example of a driver routine you might use. */
33/* Basically, the driver routine needs to call TFset_parameters() to
34set all the scalar parameters, and then call TFfit_onek() for each
35wavenumber k you desire. */
36
37/* While the routines use Mpc^-1 units internally, this driver has been
38written to take an array of wavenumbers in units of h Mpc^-1. On the
39other hand, if you want to use Mpc^-1 externally, you can do this by
40altering the variables you pass to the driver:
41 omega0 -> omega0*hubble*hubble, hubble -> 1.0 */
42
43/* INPUT: omega0 -- the matter density (baryons+CDM) in units of critical
44 f_baryon -- the ratio of baryon density to matter density
45 hubble -- the Hubble constant, in units of 100 km/s/Mpc
46 Tcmb -- the CMB temperature in Kelvin. T<=0 uses the COBE value 2.728.
47 numk -- the length of the following zero-offset array
48 k[] -- the array of wavevectors k[0..numk-1] */
49
50/* INPUT/OUTPUT: There are three output arrays of transfer functions.
51All are zero-offset and, if used, must have storage [0..numk-1] declared
52in the calling program. However, if you substitute the NULL pointer for
53one or more of the arrays, then that particular transfer function won't
54be outputted. The transfer functions are:
55
56 tf_full[] -- The full fitting formula, eq. (16), for the matter
57 transfer function.
58 tf_baryon[] -- The baryonic piece of the full fitting formula, eq. 21.
59 tf_cdm[] -- The CDM piece of the full fitting formula, eq. 17. */
60
61/* Again, you can set these pointers to NULL in the function call if
62you don't want a particular output. */
63
64/* Various intermediate scalar quantities are stored in global variables,
65so that you might more easily access them. However, this also means that
66you would be better off not simply #include'ing this file in your programs,
67but rather compiling it separately, calling only the driver, and using
68extern declarations to access the intermediate quantities. */
69
70/* ----------------------------- DRIVER ------------------------------- */
71
72void TFfit_hmpc(double omega0, double f_baryon, double hubble, double Tcmb,
73 int numk, double *k, double *tf_full, double *tf_baryon, double *tf_cdm)
74/* Remember: k[0..numk-1] is in units of h Mpc^-1. */
75{
76 int j;
77 double tf_thisk, baryon_piece, cdm_piece;
78
79 TFset_parameters(omega0*hubble*hubble, f_baryon, Tcmb);
80
81 for (j=0;j<numk;j++) {
82 tf_thisk = TFfit_onek(k[j]*hubble, &baryon_piece, &cdm_piece);
83 if (tf_full!=NULL) tf_full[j] = tf_thisk;
84 if (tf_baryon!=NULL) tf_baryon[j] = baryon_piece;
85 if (tf_cdm!=NULL) tf_cdm[j] = cdm_piece;
86 }
87 return;
88}
89
90/* ------------------------ FITTING FORMULAE ROUTINES ----------------- */
91
92/* There are two routines here. TFset_parameters() sets all the scalar
93parameters, while TFfit_onek() calculates the transfer function for a
94given wavenumber k. TFfit_onek() may be called many times after a single
95call to TFset_parameters() */
96
97/* Global variables -- We've left many of the intermediate results as
98global variables in case you wish to access them, e.g. by declaring
99them as extern variables in your main program. */
100/* Note that all internal scales are in Mpc, without any Hubble constants! */
101
102double omhh, /* Omega_matter*h^2 */
103 obhh, /* Omega_baryon*h^2 */
104 theta_cmb, /* Tcmb in units of 2.7 K */
105 z_equality, /* Redshift of matter-radiation equality, really 1+z */
106 k_equality, /* Scale of equality, in Mpc^-1 */
107 z_drag, /* Redshift of drag epoch */
108 R_drag, /* Photon-baryon ratio at drag epoch */
109 R_equality, /* Photon-baryon ratio at equality epoch */
110 sound_horizon, /* Sound horizon at drag epoch, in Mpc */
111 k_silk, /* Silk damping scale, in Mpc^-1 */
112 alpha_c, /* CDM suppression */
113 beta_c, /* CDM log shift */
114 alpha_b, /* Baryon suppression */
115 beta_b, /* Baryon envelope shift */
116 beta_node, /* Sound horizon shift */
117 k_peak, /* Fit to wavenumber of first peak, in Mpc^-1 */
118 sound_horizon_fit, /* Fit to sound horizon, in Mpc */
119 alpha_gamma; /* Gamma suppression in approximate TF */
120
121/* Convenience from Numerical Recipes in C, 2nd edition
122------------ COMPLETEMENT DEBILE (cmv) ????????? ------------
123static double sqrarg;
124#define SQR(a) ((sqrarg=(a)) == 0.0 ? 0.0 : sqrarg*sqrarg)
125static double cubearg;
126#define CUBE(a) ((cubearg=(a)) == 0.0 ? 0.0 : cubearg*cubearg*cubearg)
127static double pow4arg;
128#define POW4(a) ((pow4arg=(a)) == 0.0 ? 0.0 : pow4arg*pow4arg*pow4arg*pow4arg)
129 Yes, I know the last one isn't optimal; it doesn't appear much */
130#define SQR(a) (a*a)
131#define CUBE(a) (a*a*a)
132#define POW4(a) (a*a*a*a)
133
134void TFset_parameters(double omega0hh, double f_baryon, double Tcmb)
135/* Set all the scalars quantities for Eisenstein & Hu 1997 fitting formula */
136/* Input: omega0hh -- The density of CDM and baryons, in units of critical dens,
137 multiplied by the square of the Hubble constant, in units
138 of 100 km/s/Mpc */
139/* f_baryon -- The fraction of baryons to CDM */
140/* Tcmb -- The temperature of the CMB in Kelvin. Tcmb<=0 forces use
141 of the COBE value of 2.728 K. */
142/* Output: Nothing, but set many global variables used in TFfit_onek().
143You can access them yourself, if you want. */
144/* Note: Units are always Mpc, never h^-1 Mpc. */
145{
146 double z_drag_b1, z_drag_b2;
147 double alpha_c_a1, alpha_c_a2, beta_c_b1, beta_c_b2, alpha_b_G, y;
148
149 if (f_baryon<=0.0 || omega0hh<=0.0) {
150 fprintf(stderr, "TFset_parameters(): Illegal input.\n");
151 exit(1);
152 }
153 omhh = omega0hh;
154 obhh = omhh*f_baryon;
155 if (Tcmb<=0.0) Tcmb=2.728; /* COBE FIRAS */
156 theta_cmb = Tcmb/2.7;
157
158 z_equality = 2.50e4*omhh/POW4(theta_cmb); /* Really 1+z */
159 k_equality = 0.0746*omhh/SQR(theta_cmb);
160
161 z_drag_b1 = 0.313*pow(omhh,-0.419)*(1+0.607*pow(omhh,0.674));
162 z_drag_b2 = 0.238*pow(omhh,0.223);
163 z_drag = 1291*pow(omhh,0.251)/(1+0.659*pow(omhh,0.828))*
164 (1+z_drag_b1*pow(obhh,z_drag_b2));
165
166 R_drag = 31.5*obhh/POW4(theta_cmb)*(1000/(1+z_drag));
167 R_equality = 31.5*obhh/POW4(theta_cmb)*(1000/z_equality);
168
169 sound_horizon = 2./3./k_equality*sqrt(6./R_equality)*
170 log((sqrt(1+R_drag)+sqrt(R_drag+R_equality))/(1+sqrt(R_equality)));
171
172 k_silk = 1.6*pow(obhh,0.52)*pow(omhh,0.73)*(1+pow(10.4*omhh,-0.95));
173
174 alpha_c_a1 = pow(46.9*omhh,0.670)*(1+pow(32.1*omhh,-0.532));
175 alpha_c_a2 = pow(12.0*omhh,0.424)*(1+pow(45.0*omhh,-0.582));
176 alpha_c = pow(alpha_c_a1,-f_baryon)*
177 pow(alpha_c_a2,-CUBE(f_baryon));
178
179 beta_c_b1 = 0.944/(1+pow(458*omhh,-0.708));
180 beta_c_b2 = pow(0.395*omhh, -0.0266);
181 beta_c = 1.0/(1+beta_c_b1*(pow(1-f_baryon, beta_c_b2)-1));
182
183 y = z_equality/(1+z_drag);
184 alpha_b_G = y*(-6.*sqrt(1+y)+(2.+3.*y)*log((sqrt(1+y)+1)/(sqrt(1+y)-1)));
185 alpha_b = 2.07*k_equality*sound_horizon*pow(1+R_drag,-0.75)*alpha_b_G;
186
187 beta_node = 8.41*pow(omhh, 0.435);
188 beta_b = 0.5+f_baryon+(3.-2.*f_baryon)*sqrt(pow(17.2*omhh,2.0)+1);
189
190 k_peak = 2.5*3.14159*(1+0.217*omhh)/sound_horizon;
191 sound_horizon_fit = 44.5*log(9.83/omhh)/sqrt(1+10.0*pow(obhh,0.75));
192
193 alpha_gamma = 1-0.328*log(431.0*omhh)*f_baryon + 0.38*log(22.3*omhh)*
194 SQR(f_baryon);
195
196 return;
197}
198
199double TFfit_onek(double k, double *tf_baryon, double *tf_cdm)
200/* Input: k -- Wavenumber at which to calculate transfer function, in Mpc^-1.
201 *tf_baryon, *tf_cdm -- Input value not used; replaced on output if
202 the input was not NULL. */
203/* Output: Returns the value of the full transfer function fitting formula.
204 This is the form given in Section 3 of Eisenstein & Hu (1997).
205 *tf_baryon -- The baryonic contribution to the full fit.
206 *tf_cdm -- The CDM contribution to the full fit. */
207/* Notes: Units are Mpc, not h^-1 Mpc. */
208{
209 double T_c_ln_beta, T_c_ln_nobeta, T_c_C_alpha, T_c_C_noalpha;
210 double q, xx, xx_tilde, q_eff;
211 double T_c_f, T_c, s_tilde, T_b_T0, T_b, f_baryon, T_full;
212 double T_0_L0, T_0_C0, T_0, gamma_eff;
213 double T_nowiggles_L0, T_nowiggles_C0, T_nowiggles;
214
215 k = fabs(k); /* Just define negative k as positive */
216 if (k==0.0) {
217 if (tf_baryon!=NULL) *tf_baryon = 1.0;
218 if (tf_cdm!=NULL) *tf_cdm = 1.0;
219 return 1.0;
220 }
221
222 q = k/13.41/k_equality;
223 xx = k*sound_horizon;
224
225 T_c_ln_beta = log(2.718282+1.8*beta_c*q);
226 T_c_ln_nobeta = log(2.718282+1.8*q);
227 T_c_C_alpha = 14.2/alpha_c + 386.0/(1+69.9*pow(q,1.08));
228 T_c_C_noalpha = 14.2 + 386.0/(1+69.9*pow(q,1.08));
229
230 T_c_f = 1.0/(1.0+POW4(xx/5.4));
231 T_c = T_c_f*T_c_ln_beta/(T_c_ln_beta+T_c_C_noalpha*SQR(q)) +
232 (1-T_c_f)*T_c_ln_beta/(T_c_ln_beta+T_c_C_alpha*SQR(q));
233
234 s_tilde = sound_horizon*pow(1+CUBE(beta_node/xx),-1./3.);
235 xx_tilde = k*s_tilde;
236
237 T_b_T0 = T_c_ln_nobeta/(T_c_ln_nobeta+T_c_C_noalpha*SQR(q));
238 T_b = sin(xx_tilde)/(xx_tilde)*(T_b_T0/(1+SQR(xx/5.2))+
239 alpha_b/(1+CUBE(beta_b/xx))*exp(-pow(k/k_silk,1.4)));
240
241 f_baryon = obhh/omhh;
242 T_full = f_baryon*T_b + (1-f_baryon)*T_c;
243
244 /* Now to store these transfer functions */
245 if (tf_baryon!=NULL) *tf_baryon = T_b;
246 if (tf_cdm!=NULL) *tf_cdm = T_c;
247 return T_full;
248}
249
250/* ======================= Approximate forms =========================== */
251
252double TFsound_horizon_fit(double omega0, double f_baryon, double hubble)
253/* Input: omega0 -- CDM density, in units of critical density
254 f_baryon -- Baryon fraction, the ratio of baryon to CDM density.
255 hubble -- Hubble constant, in units of 100 km/s/Mpc */
256/* Output: The approximate value of the sound horizon, in h^-1 Mpc. */
257/* Note: If you prefer to have the answer in units of Mpc, use hubble -> 1
258and omega0 -> omega0*hubble^2. */
259{
260 double omhh, sound_horizon_fit_mpc;
261 omhh = omega0*hubble*hubble;
262 sound_horizon_fit_mpc =
263 44.5*log(9.83/omhh)/sqrt(1+10.0*pow(omhh*f_baryon,0.75));
264 return sound_horizon_fit_mpc*hubble;
265}
266
267double TFk_peak(double omega0, double f_baryon, double hubble)
268/* Input: omega0 -- CDM density, in units of critical density
269 f_baryon -- Baryon fraction, the ratio of baryon to CDM density.
270 hubble -- Hubble constant, in units of 100 km/s/Mpc */
271/* Output: The approximate location of the first baryonic peak, in h Mpc^-1 */
272/* Note: If you prefer to have the answer in units of Mpc^-1, use hubble -> 1
273and omega0 -> omega0*hubble^2. */
274{
275 double omhh, k_peak_mpc;
276 omhh = omega0*hubble*hubble;
277 k_peak_mpc = 2.5*3.14159*(1+0.217*omhh)/TFsound_horizon_fit(omhh,f_baryon,1.0);
278 return k_peak_mpc/hubble;
279}
280
281double TFnowiggles(double omega0, double f_baryon, double hubble,
282 double Tcmb, double k_hmpc)
283/* Input: omega0 -- CDM density, in units of critical density
284 f_baryon -- Baryon fraction, the ratio of baryon to CDM density.
285 hubble -- Hubble constant, in units of 100 km/s/Mpc
286 Tcmb -- Temperature of the CMB in Kelvin; Tcmb<=0 forces use of
287 COBE FIRAS value of 2.728 K
288 k_hmpc -- Wavenumber in units of (h Mpc^-1). */
289/* Output: The value of an approximate transfer function that captures the
290non-oscillatory part of a partial baryon transfer function. In other words,
291the baryon oscillations are left out, but the suppression of power below
292the sound horizon is included. See equations (30) and (31). */
293/* Note: If you prefer to use wavenumbers in units of Mpc^-1, use hubble -> 1
294and omega0 -> omega0*hubble^2. */
295{
296 double k, omhh, theta_cmb, k_equality, q, xx, alpha_gamma, gamma_eff;
297 double q_eff, T_nowiggles_L0, T_nowiggles_C0;
298
299 k = k_hmpc*hubble; /* Convert to Mpc^-1 */
300 omhh = omega0*hubble*hubble;
301 if (Tcmb<=0.0) Tcmb=2.728; /* COBE FIRAS */
302 theta_cmb = Tcmb/2.7;
303
304 k_equality = 0.0746*omhh/SQR(theta_cmb);
305 q = k/13.41/k_equality;
306 xx = k*TFsound_horizon_fit(omhh, f_baryon, 1.0);
307
308 alpha_gamma = 1-0.328*log(431.0*omhh)*f_baryon + 0.38*log(22.3*omhh)*
309 SQR(f_baryon);
310 gamma_eff = omhh*(alpha_gamma+(1-alpha_gamma)/(1+POW4(0.43*xx)));
311 q_eff = q*omhh/gamma_eff;
312
313 T_nowiggles_L0 = log(2.0*2.718282+1.8*q_eff);
314 T_nowiggles_C0 = 14.2 + 731.0/(1+62.5*q_eff);
315 return T_nowiggles_L0/(T_nowiggles_L0+T_nowiggles_C0*SQR(q_eff));
316}
317
318/* ======================= Zero Baryon Formula =========================== */
319
320double TFzerobaryon(double omega0, double hubble, double Tcmb, double k_hmpc)
321/* Input: omega0 -- CDM density, in units of critical density
322 hubble -- Hubble constant, in units of 100 km/s/Mpc
323 Tcmb -- Temperature of the CMB in Kelvin; Tcmb<=0 forces use of
324 COBE FIRAS value of 2.728 K
325 k_hmpc -- Wavenumber in units of (h Mpc^-1). */
326/* Output: The value of the transfer function for a zero-baryon universe. */
327/* Note: If you prefer to use wavenumbers in units of Mpc^-1, use hubble -> 1
328and omega0 -> omega0*hubble^2. */
329{
330 double k, omhh, theta_cmb, k_equality, q, T_0_L0, T_0_C0;
331
332 k = k_hmpc*hubble; /* Convert to Mpc^-1 */
333 omhh = omega0*hubble*hubble;
334 if (Tcmb<=0.0) Tcmb=2.728; /* COBE FIRAS */
335 theta_cmb = Tcmb/2.7;
336
337 k_equality = 0.0746*omhh/SQR(theta_cmb);
338 q = k/13.41/k_equality;
339
340 T_0_L0 = log(2.0*2.718282+1.8*q);
341 T_0_C0 = 14.2 + 731.0/(1+62.5*q);
342 return T_0_L0/(T_0_L0+T_0_C0*q*q);
343}
344
345
346/* Added by cmv 24/08/2007 */
347void TFprint_parameters(void)
348{
349 printf("omhh=%g obhh=%g theta_cmb=%g\n",omhh,obhh,theta_cmb);
350 printf("z_equality=%g k_equality=%g R_equality=%g\n",z_equality,k_equality,R_equality);
351 printf("z_drag=%g R_drag=%g sound_horizon=%g\n",z_drag,R_drag,sound_horizon);
352 printf("k_silk=%g\n",k_silk);
353 printf("alpha_c=%g beta_c=%g\n",alpha_c,beta_c);
354 printf("alpha_b=%g beta_b=%g\n",alpha_b,beta_b);
355 printf("beta_node=%g\n",beta_node);
356 printf("k_peak=%g sound_horizon_fit=%g\n",k_peak,sound_horizon_fit);
357 printf("alpha_gamma=%g\n",alpha_gamma);
358}
Note: See TracBrowser for help on using the repository browser.