00001 #include <math.h>
00002 #include <gsl/gsl_integration.h>
00003 #include <gsl/gsl_sf_gamma.h>
00004 #include "Cosmology.h"
00005
00006 const double D_H_w_h_Mpc = 2997.925;
00007 const double t_H_w_h_year = 9.78e+9;
00008
00009 Cosmology::Cosmology(double o_m, double o_x, double o_k, \
00010 double eos_w_x, double hubble)
00011 {
00012 omega_m = o_m;
00013 omega_x = o_x;
00014 omega_k = o_k;
00015 w_x = eos_w_x;
00016 h = hubble;
00017 }
00018
00019 Cosmology::~Cosmology()
00020 {
00021 }
00022
00023 struct E_z_par {
00024 double omega_m;
00025 double omega_x;
00026 double omega_k;
00027 double w_x;
00028 double h;
00029 };
00030
00031 double E_z(double x, void *params) {
00032 double omega_m, omega_x, omega_k, w_x, h;
00033 double fun;
00034 struct E_z_par *temp;
00035 temp = (struct E_z_par *) params;
00036
00037 omega_m = temp->omega_m;
00038 omega_x = temp->omega_x;
00039 omega_k = temp->omega_k;
00040 w_x = temp->w_x;
00041 h = temp->h;
00042
00043 if(w_x == -1.0) {
00044 fun = sqrt(omega_m*(1.0 + x)*(1.0 + x)*(1.0 + x) + \
00045 omega_k*(1.0 + x)*(1.0 + x) + omega_x);
00046 } else {
00047 fun = sqrt(omega_m*(1.0 + x)*(1.0 + x)*(1.0 + x) + \
00048 omega_k*(1.0 + x)*(1.0 + x) + \
00049 omega_x*pow((1.0 + x), 3.0*(1.0 + w_x)));
00050 }
00051 fun = 1.0/fun;
00052 return fun;
00053 }
00054
00055 double Cosmology::D_C(double z1, double z2)
00056 {
00057 double return_val;
00058 gsl_integration_workspace *w = gsl_integration_workspace_alloc(2000);
00059 double result, error;
00060 gsl_function F;
00061 struct E_z_par temp;
00062
00063 temp.omega_m = omega_m;
00064 temp.omega_x = omega_x;
00065 temp.omega_k = omega_k;
00066 temp.w_x = w_x;
00067 temp.h = h;
00068
00069 F.function = &E_z;
00070 F.params = &temp;
00071 gsl_integration_qag(&F, z1, z2, 1.0e-7, 1.0e-8, 2000, 6, w, &result, &error);
00072 return_val = (D_H_w_h_Mpc / h) * result;
00073 gsl_integration_workspace_free(w);
00074
00075 return return_val;
00076 }
00077
00078 double Cosmology::D_M(double z1, double z2)
00079 {
00080 double return_val;
00081
00082 if( omega_k == 0.0) {
00083 return_val = D_C(z1, z2);
00084 return return_val;
00085 }
00086
00087 if( omega_k > 0.0) {
00088 return_val = sinh(sqrt(omega_k) * D_C(z1, z2) * h / D_H_w_h_Mpc);
00089 return_val = D_H_w_h_Mpc * return_val / (h * sqrt(omega_k));
00090 return return_val;
00091 }
00092
00093 if( omega_k < 0.0) {
00094 return_val = sin(sqrt(omega_k) * D_C(z1, z2) * h / D_H_w_h_Mpc);
00095 return_val = D_H_w_h_Mpc * return_val / (h * sqrt(omega_k));
00096 return return_val;
00097 }
00098
00099 return -1.0;
00100 }
00101
00102 double Cosmology::ang_dist(double z)
00103 {
00104 double return_val;
00105
00106 return_val = D_M(0.0, z)/(1.0 + z);
00107 return return_val;
00108 }
00109
00110 double Cosmology::ang_dist_z1_z2(double z1, double z2)
00111 {
00112 double return_val;
00113
00114 return_val = D_M(z1, z2) / (1.0 + z2);
00115 return return_val;
00116 }
00117
00118 double Cosmology::lum_dist(double z)
00119 {
00120 double return_val;
00121
00122 return_val = D_M(0.0, z) * (1.0 + z);
00123 return return_val;
00124 }
00125
00126 double Cosmology::unit_comoving_volume(double z, void *params)
00127 {
00128 double D_A;
00129 double return_val;
00130 struct E_z_par temp;
00131
00132 temp.omega_m = omega_m;
00133 temp.omega_x = omega_x;
00134 temp.omega_k = omega_k;
00135 temp.w_x = w_x;
00136 temp.h = h;
00137
00138 D_A = ang_dist(z);
00139 return_val = (D_H_w_h_Mpc / h) * (1.0 + z) * (1.0 + z);
00140 return_val = return_val * D_A * D_A * E_z(z, &temp);
00141
00142 return return_val;
00143 }
00144
00145 double comoving_volume_integral(double z, void *params)
00146 {
00147 double D_A;
00148 double return_val;
00149 struct E_z_par *temp_pt;
00150 struct E_z_par temp;
00151 Cosmology tmp;
00152 temp_pt = (struct E_z_par *) params;
00153
00154 temp.omega_m = temp_pt->omega_m;
00155 temp.omega_x = temp_pt->omega_x;
00156 temp.omega_k = temp_pt->omega_k;
00157 temp.w_x = temp_pt->w_x;
00158 temp.h = temp_pt->h;
00159
00160 tmp = Cosmology(temp.omega_m, temp.omega_x, \
00161 temp.omega_k, temp.w_x, temp.h);
00162 D_A = tmp.ang_dist(z);
00163 return_val = (D_H_w_h_Mpc / temp.h) * (1.0 + z) * (1.0 + z);
00164 return_val = return_val * D_A * D_A * E_z(z, &temp);
00165
00166 return return_val;
00167 }
00168
00169
00170 double Cosmology::comoving_volume_z1_z2(double z1, double z2)
00171 {
00172 gsl_integration_workspace *w = gsl_integration_workspace_alloc(2000);
00173 double result, error;
00174 gsl_function F;
00175 struct E_z_par temp;
00176
00177 temp.omega_m = omega_m;
00178 temp.omega_x = omega_x;
00179 temp.omega_k = omega_k;
00180 temp.w_x = w_x;
00181 temp.h = h;
00182
00183 F.function = &comoving_volume_integral;
00184 F.params = &temp;
00185 gsl_integration_qag(&F, z1, z2, 1.0e-7, 1.0e-8, 2000, 6, w, &result, &error);
00186 gsl_integration_workspace_free(w);
00187
00188 return result;
00189 }
00190
00191 double Cosmology::comoving_volume_0_z(double z)
00192 {
00193 double return_val;
00194
00195 if(omega_k > 0.0) {
00196 return_val = (D_M(0.0, z) * h / D_H_w_h_Mpc) * \
00197 sqrt(1.0 + omega_k * (D_M(0.0, z) * h / D_H_w_h_Mpc) * \
00198 (D_M(0.0, z) * h / D_H_w_h_Mpc));
00199 return_val = return_val - (1.0/sqrt(fabs(omega_k))) * \
00200 asinh(sqrt(fabs(omega_k)) * (D_M(0.0, z) * h / D_H_w_h_Mpc));
00201 return_val = (return_val / 2.0) * \
00202 (D_H_w_h_Mpc * D_H_w_h_Mpc * D_H_w_h_Mpc) / omega_k;
00203 return_val = return_val / (h * h * h);
00204 return return_val;
00205 }
00206
00207 if(omega_k < 0.0) {
00208 return_val = (D_M(0.0, z) * h / D_H_w_h_Mpc) * \
00209 sqrt(1.0 + omega_k * (D_M(0.0, z) * h / D_H_w_h_Mpc) * \
00210 (D_M(0.0, z) * h / D_H_w_h_Mpc));
00211 return_val = return_val - (1.0/sqrt(fabs(omega_k))) * \
00212 asin(sqrt(fabs(omega_k)) * (D_M(0.0, z) * h / D_H_w_h_Mpc));
00213 return_val = (return_val / 2.0) * \
00214 (D_H_w_h_Mpc * D_H_w_h_Mpc * D_H_w_h_Mpc) / omega_k;
00215 return_val = return_val / (h * h * h);
00216 return return_val;
00217 }
00218
00219
00220 return_val = D_M(0.0, z) * D_M(0.0, z) * D_M(0.0, z) / 3.0;
00221
00222 return return_val;
00223 }
00224
00225 double time_integral(double x, void *params)
00226 {
00227 double fun;
00228 struct E_z_par* temp_pt;
00229 struct E_z_par temp;
00230 temp_pt = (struct E_z_par *) params;
00231
00232 temp.omega_m = temp_pt->omega_m;
00233 temp.omega_x = temp_pt->omega_x;
00234 temp.omega_k = temp_pt->omega_k;
00235 temp.w_x = temp_pt->w_x;
00236 temp.h = temp_pt->h;
00237
00238 fun = (1.0 / (1.0 + x)) * E_z(x, &temp);
00239
00240 return fun;
00241 }
00242
00243 double Cosmology::lookback_time(double z)
00244 {
00245 gsl_integration_workspace *w = gsl_integration_workspace_alloc(2000);
00246 double result, error;
00247 gsl_function F;
00248 struct E_z_par temp;
00249
00250 temp.omega_m = omega_m;
00251 temp.omega_x = omega_x;
00252 temp.omega_k = omega_k;
00253 temp.w_x = w_x;
00254 temp.h = h;
00255
00256 F.function = &time_integral;
00257 F.params = &temp;
00258 gsl_integration_qag(&F, 0.0, z, 1.0e-7, 1.0e-8, 2000, 6, w, &result, &error);
00259 result = result * t_H_w_h_year / h;
00260 gsl_integration_workspace_free(w);
00261
00262 return result;
00263 }
00264
00265 double Cosmology::age(double z)
00266 {
00267 gsl_integration_workspace *w = gsl_integration_workspace_alloc(2000);
00268 double result, error;
00269 gsl_function F;
00270 struct E_z_par temp;
00271
00272 temp.omega_m = omega_m;
00273 temp.omega_x = omega_x;
00274 temp.omega_k = omega_k;
00275 temp.w_x = w_x;
00276 temp.h = h;
00277
00278 F.function = &time_integral;
00279 F.params = &temp;
00280 gsl_integration_qagiu(&F, z, 1.0e-7, 1.0e-8, 2000, w, &result, &error);
00281 result = result * t_H_w_h_year / h;
00282 gsl_integration_workspace_free(w);
00283
00284 return result;
00285 }
00286
00287 double Cosmology::age_now(void)
00288 {
00289 double return_val;
00290
00291 return_val = age(0.0);
00292
00293 return return_val;
00294 }
00295
00296 double schechter(double L, double phi_star, double alpha, double L_star)
00297 {
00298 double return_val;
00299 L_star = 1.0 / L_star;
00300
00301 return_val = phi_star;
00302 return_val = return_val * pow(L*L_star, alpha);
00303 return_val = return_val * exp(-1.0*L*L_star) * L_star;
00304
00305 return return_val;
00306 }
00307
00308 double number_schechter(double L, double phi_star, double alpha, double L_star)
00309 {
00310 double return_val;
00311
00312 return_val = phi_star * gsl_sf_gamma_inc(alpha+1.0, L/L_star);
00313
00314 return return_val;
00315 }
00316
00317 double total_number_schechter(double phi_star, double alpha)
00318 {
00319 double return_val;
00320
00321 return_val = phi_star * gsl_sf_gamma(alpha+1.0);
00322
00323 return return_val;
00324 }
00325
00326 double luminosity_schechter(double L, double phi_star, double alpha, double L_star)
00327 {
00328 double return_val;
00329
00330 return_val = phi_star * L_star * gsl_sf_gamma_inc(alpha+2.0, L/L_star);
00331
00332 return return_val;
00333 }
00334
00335 double total_luminosity_schechter(double phi_star, double alpha, double L_star)
00336 {
00337 double return_val;
00338
00339 return_val = phi_star * L_star * gsl_sf_gamma(alpha+2.0);
00340
00341 return return_val;
00342 }
00343
00344 double double_schechter(double L, double phi_star1, double alpha1, \
00345 double L_star1, double phi_star2, double alpha2, double L_star2)
00346 {
00347 double return_val, temp;
00348
00349 L_star1 = 1.0 / L_star1;
00350 return_val = phi_star1;
00351 return_val = return_val * pow(L*L_star1, alpha1);
00352 return_val = return_val * exp(-1.0*L*L_star1) * L_star1;
00353 L_star2 = 1.0 / L_star2;
00354 temp = phi_star2;
00355 temp = temp * pow(L*L_star2, alpha2);
00356 temp = temp * exp(-1.0*L*L_star2) * L_star2;
00357
00358 return_val = return_val + temp;
00359
00360 return return_val;
00361 }
00362
00363 double number_double_schechter(double L, double phi_star1, double alpha1, double L_star1, \
00364 double phi_star2, double alpha2, double L_star2)
00365 {
00366 double return_val;
00367
00368 return_val = phi_star1 * gsl_sf_gamma_inc(alpha1+1.0, L/L_star1);
00369 return_val = return_val + phi_star2 * gsl_sf_gamma_inc(alpha2+1.0, L/L_star2);
00370
00371 return return_val;
00372 }
00373
00374 double luminosity_double_schechter(double L, double phi_star1, double alpha1, double L_star1, \
00375 double phi_star2, double alpha2, double L_star2)
00376 {
00377 double return_val;
00378
00379 return_val = phi_star1 * L_star1 * gsl_sf_gamma_inc(alpha1+2.0, L/L_star1);
00380 return_val = return_val + phi_star2 * L_star2 * gsl_sf_gamma_inc(alpha2+2.0, L/L_star2);
00381
00382 return return_val;
00383 }
00384
00385 double total_luminosity_double_schechter(double phi_star1, double alpha1, double L_star1, \
00386 double phi_star2, double alpha2, double L_star2)
00387 {
00388 double return_val;
00389
00390 return_val = phi_star1 * L_star1 * gsl_sf_gamma(alpha1+2.0);
00391 return_val = return_val + phi_star2 * L_star2 * gsl_sf_gamma(alpha2+2.0);
00392
00393 return return_val;
00394 }