Package astLib :: Module astCalc
[hide private]
[frames] | no frames]

Source Code for Module astLib.astCalc

  1  """module for performing common calculations 
  2   
  3  (c) 2007-2011 Matt Hilton 
  4   
  5  (c) 2013-2014 Matt Hilton & Steven Boada 
  6   
  7  U{http://astlib.sourceforge.net} 
  8   
  9  The focus in this module is at present on calculations of distances in a given 
 10  cosmology. The parameters for the cosmological model are set using the 
 11  variables OMEGA_M0, OMEGA_L0, OMEGA_R0, H0 in the module namespace (see below for details). 
 12   
 13  @var OMEGA_M0: The matter density parameter at z=0. 
 14  @type OMEGA_M0: float 
 15   
 16  @var OMEGA_L0: The dark energy density (in the form of a cosmological 
 17      constant) at z=0. 
 18  @type OMEGA_L0: float 
 19   
 20  @var OMEGA_R0: The radiation density at z=0 (note this is only used currently 
 21      in calculation of L{Ez}). 
 22  @type OMEGA_R0: float 
 23   
 24  @var H0: The Hubble parameter (in km/s/Mpc) at z=0. 
 25  @type H0: float 
 26   
 27  @var C_LIGHT: The speed of light in km/s. 
 28  @type C_LIGHT: float 
 29   
 30  """ 
 31   
 32  OMEGA_M0 = 0.3 
 33  OMEGA_L0 = 0.7 
 34  OMEGA_R0 = 8.24E-5 
 35  H0 = 70.0 
 36   
 37  C_LIGHT = 3.0e5 
 38   
 39  import math 
 40  try: 
 41      from scipy import integrate 
 42  except ImportError: 
 43      print "WARNING: astCalc failed to import scipy modules - ", 
 44      print "some functions will not work" 
 45   
 46  #import sys 
 47   
 48  #------------------------------------------------------------------------------ 
49 -def dl(z):
50 """Calculates the luminosity distance in Mpc at redshift z. 51 52 @type z: float 53 @param z: redshift 54 @rtype: float 55 @return: luminosity distance in Mpc 56 57 """ 58 59 DM = dm(z) 60 DL = (1.0+z)*DM 61 62 return DL
63 64 #------------------------------------------------------------------------------
65 -def da(z):
66 """Calculates the angular diameter distance in Mpc at redshift z. 67 68 @type z: float 69 @param z: redshift 70 @rtype: float 71 @return: angular diameter distance in Mpc 72 73 """ 74 DM = dm(z) 75 DA = DM/(1.0+z) 76 77 return DA
78 79 #------------------------------------------------------------------------------
80 -def dm(z):
81 """Calculates the transverse comoving distance (proper motion distance) in 82 Mpc at redshift z. 83 84 @type z: float 85 @param z: redshift 86 @rtype: float 87 @return: transverse comoving distance (proper motion distance) in Mpc 88 89 """ 90 91 OMEGA_K = 1.0 - OMEGA_M0 - OMEGA_L0 92 93 # Integration limits 94 xMax = 1.0 95 xMin = 1.0 / (1.0 + z) 96 97 # Function to be integrated 98 yn = lambda x: (1.0/math.sqrt(OMEGA_M0*x + OMEGA_L0*math.pow(x, 4) + 99 OMEGA_K*math.pow(x, 2))) 100 101 integralValue, integralError = integrate.quad(yn, xMin, xMax) 102 103 if OMEGA_K > 0.0: 104 DM = (C_LIGHT/H0 * math.pow(abs(OMEGA_K), -0.5) * 105 math.sinh(math.sqrt(abs(OMEGA_K)) * integralValue)) 106 elif OMEGA_K == 0.0: 107 DM = C_LIGHT/H0 * integralValue 108 elif OMEGA_K < 0.0: 109 DM = (C_LIGHT/H0 * math.pow(abs(OMEGA_K), -0.5) * 110 math.sin(math.sqrt(abs(OMEGA_K)) * integralValue)) 111 112 return DM
113 114 #------------------------------------------------------------------------------
115 -def dc(z):
116 """Calculates the line of sight comoving distance in Mpc at redshift z. 117 118 @type z: float 119 @param z: redshift 120 @rtype: float 121 @return: transverse comoving distance (proper motion distance) in Mpc 122 123 """ 124 125 OMEGA_K = 1.0 - OMEGA_M0 - OMEGA_L0 126 127 # Integration limits 128 xMax = 1.0 129 xMin = 1.0 / (1.0 + z) 130 131 # Function to be integrated 132 yn = lambda x: (1.0/math.sqrt(OMEGA_M0*x + OMEGA_L0*math.pow(x, 4) + 133 OMEGA_K*math.pow(x, 2))) 134 135 integralValue, integralError = integrate.quad(yn, xMin, xMax) 136 137 DC= C_LIGHT/H0*integralValue 138 139 return DC
140 141 #------------------------------------------------------------------------------
142 -def dVcdz(z):
143 """Calculates the line of sight comoving volume element per steradian dV/dz 144 at redshift z. 145 146 @type z: float 147 @param z: redshift 148 @rtype: float 149 @return: comoving volume element per steradian 150 151 """ 152 153 dH = C_LIGHT/H0 154 dVcdz=(dH*(math.pow(da(z),2))*(math.pow(1+z,2))/Ez(z)) 155 156 return dVcdz
157 158 #------------------------------------------------------------------------------
159 -def dl2z(distanceMpc):
160 """Calculates the redshift z corresponding to the luminosity distance given 161 in Mpc. 162 163 @type distanceMpc: float 164 @param distanceMpc: distance in Mpc 165 @rtype: float 166 @return: redshift 167 168 """ 169 170 dTarget = distanceMpc 171 172 toleranceMpc = 0.1 173 174 zMin = 0.0 175 zMax = 10.0 176 177 diff = dl(zMax) - dTarget 178 while diff < 0: 179 zMax = zMax + 5.0 180 diff = dl(zMax) - dTarget 181 182 zTrial = zMin + (zMax-zMin)/2.0 183 184 dTrial = dl(zTrial) 185 diff = dTrial - dTarget 186 while abs(diff) > toleranceMpc: 187 188 if diff > 0: 189 zMax = zMax - (zMax-zMin)/2.0 190 else: 191 zMin = zMin + (zMax-zMin)/2.0 192 193 zTrial = zMin + (zMax-zMin)/2.0 194 dTrial = dl(zTrial) 195 diff = dTrial - dTarget 196 197 return zTrial
198 199 #------------------------------------------------------------------------------
200 -def dc2z(distanceMpc):
201 """Calculates the redshift z corresponding to the comoving distance given 202 in Mpc. 203 204 @type distanceMpc: float 205 @param distanceMpc: distance in Mpc 206 @rtype: float 207 @return: redshift 208 209 """ 210 211 dTarget = distanceMpc 212 213 toleranceMpc = 0.1 214 215 zMin = 0.0 216 zMax = 10.0 217 218 diff = dc(zMax) - dTarget 219 while diff < 0: 220 zMax = zMax + 5.0 221 diff = dc(zMax) - dTarget 222 223 zTrial = zMin + (zMax-zMin)/2.0 224 225 dTrial = dc(zTrial) 226 diff = dTrial - dTarget 227 while abs(diff) > toleranceMpc: 228 229 if diff > 0: 230 zMax = zMax - (zMax-zMin)/2.0 231 else: 232 zMin = zMin + (zMax-zMin)/2.0 233 234 zTrial = zMin + (zMax-zMin)/2.0 235 dTrial = dc(zTrial) 236 diff = dTrial - dTarget 237 238 return zTrial
239 240 #------------------------------------------------------------------------------
241 -def t0():
242 """Calculates the age of the universe in Gyr at z=0 for the current set of 243 cosmological parameters. 244 245 @rtype: float 246 @return: age of the universe in Gyr at z=0 247 248 """ 249 250 OMEGA_K = 1.0 - OMEGA_M0 - OMEGA_L0 251 252 # Integration limits 253 xMax = 1.0 254 xMin = 0 255 256 # Function to be integrated 257 yn = lambda x: (x/math.sqrt(OMEGA_M0*x + OMEGA_L0*math.pow(x, 4) + 258 OMEGA_K*math.pow(x, 2))) 259 260 integralValue, integralError = integrate.quad(yn, xMin, xMax) 261 262 T0 = (1.0/H0*integralValue*3.08e19)/3.16e7/1e9 263 264 return T0
265 266 #------------------------------------------------------------------------------
267 -def tl(z):
268 """ Calculates the lookback time in Gyr to redshift z for the current set 269 of cosmological parameters. 270 271 @type z: float 272 @param z: redshift 273 @rtype: float 274 @return: lookback time in Gyr to redshift z 275 276 """ 277 OMEGA_K = 1.0 - OMEGA_M0 - OMEGA_L0 278 279 # Integration limits 280 xMax = 1.0 281 xMin = 1./(1.+z) 282 283 # Function to be integrated 284 yn = lambda x: (x/math.sqrt(OMEGA_M0*x + OMEGA_L0*math.pow(x, 4) + 285 OMEGA_K*math.pow(x, 2))) 286 287 integralValue, integralError = integrate.quad(yn, xMin, xMax) 288 289 T0 = (1.0/H0*integralValue*3.08e19)/3.16e7/1e9 290 291 return T0
292 293 #------------------------------------------------------------------------------
294 -def tz(z):
295 """Calculates the age of the universe at redshift z for the current set of 296 cosmological parameters. 297 298 @type z: float 299 @param z: redshift 300 @rtype: float 301 @return: age of the universe in Gyr at redshift z 302 303 """ 304 305 TZ = t0() - tl(z) 306 307 return TZ
308 309 #------------------------------------------------------------------------------
310 -def tl2z(tlGyr):
311 """Calculates the redshift z corresponding to lookback time tlGyr given in 312 Gyr. 313 314 @type tlGyr: float 315 @param tlGyr: lookback time in Gyr 316 @rtype: float 317 @return: redshift 318 319 @note: Raises ValueError if tlGyr is not positive. 320 321 """ 322 if tlGyr < 0.: 323 raise ValueError('Lookback time must be positive') 324 325 tTarget = tlGyr 326 327 toleranceGyr = 0.001 328 329 zMin = 0.0 330 zMax = 10.0 331 332 diff = tl(zMax) - tTarget 333 while diff < 0: 334 zMax = zMax + 5.0 335 diff = tl(zMax) - tTarget 336 337 zTrial = zMin + (zMax-zMin)/2.0 338 339 tTrial = tl(zTrial) 340 diff = tTrial - tTarget 341 while abs(diff) > toleranceGyr: 342 343 if diff > 0: 344 zMax = zMax - (zMax-zMin)/2.0 345 else: 346 zMin = zMin + (zMax-zMin)/2.0 347 348 zTrial = zMin + (zMax-zMin)/2.0 349 tTrial = tl(zTrial) 350 diff = tTrial - tTarget 351 352 return zTrial
353 354 #------------------------------------------------------------------------------
355 -def tz2z(tzGyr):
356 """Calculates the redshift z corresponding to age of the universe tzGyr 357 given in Gyr. 358 359 @type tzGyr: float 360 @param tzGyr: age of the universe in Gyr 361 @rtype: float 362 @return: redshift 363 364 @note: Raises ValueError if Universe age not positive 365 366 """ 367 if tzGyr <= 0: 368 raise ValueError('Universe age must be positive.') 369 tl = t0() - tzGyr 370 z = tl2z(tl) 371 372 return z
373 374 #------------------------------------------------------------------------------
375 -def absMag(appMag, distMpc):
376 """Calculates the absolute magnitude of an object at given luminosity 377 distance in Mpc. 378 379 @type appMag: float 380 @param appMag: apparent magnitude of object 381 @type distMpc: float 382 @param distMpc: distance to object in Mpc 383 @rtype: float 384 @return: absolute magnitude of object 385 386 """ 387 absMag = appMag - (5.0*math.log10(distMpc*1.0e5)) 388 389 return absMag
390 391 #------------------------------------------------------------------------------
392 -def Ez(z):
393 """Calculates the value of E(z), which describes evolution of the Hubble 394 parameter with redshift, at redshift z for the current set of cosmological 395 parameters. See, e.g., Bryan & Norman 1998 (ApJ, 495, 80). 396 397 @type z: float 398 @param z: redshift 399 @rtype: float 400 @return: value of E(z) at redshift z 401 402 """ 403 404 Ez = math.sqrt(Ez2(z)) 405 406 return Ez
407 408 #------------------------------------------------------------------------------
409 -def Ez2(z):
410 """Calculates the value of E(z)^2, which describes evolution of the Hubble 411 parameter with redshift, at redshift z for the current set of cosmological 412 parameters. See, e.g., Bryan & Norman 1998 (ApJ, 495, 80). 413 414 @type z: float 415 @param z: redshift 416 @rtype: float 417 @return: value of E(z)^2 at redshift z 418 419 """ 420 # This form of E(z) is more reliable at high redshift. It is basically the 421 # same for all redshifts below 10. But above that, the radiation term 422 # begins to dominate. From Peebles 1993. 423 424 Ez2 = (OMEGA_R0 * math.pow(1.0+z, 4) + 425 OMEGA_M0* math.pow(1.0+z, 3) + 426 (1.0- OMEGA_M0- OMEGA_L0) * 427 math.pow(1.0+z, 2) + OMEGA_L0) 428 429 return Ez2
430 431 #------------------------------------------------------------------------------
432 -def OmegaMz(z):
433 """Calculates the matter density of the universe at redshift z. See, e.g., 434 Bryan & Norman 1998 (ApJ, 495, 80). 435 436 @type z: float 437 @param z: redshift 438 @rtype: float 439 @return: matter density of universe at redshift z 440 441 """ 442 ez2 = Ez2(z) 443 444 Omega_Mz = (OMEGA_M0*math.pow(1.0+z, 3))/ez2 445 446 return Omega_Mz
447 448 #------------------------------------------------------------------------------
449 -def OmegaLz(z):
450 """ Calculates the dark energy density of the universe at redshift z. 451 452 @type z: float 453 @param z: redshift 454 @rtype: float 455 @return: dark energy density of universe at redshift z 456 457 """ 458 ez2 = Ez2(z) 459 460 return OMEGA_L0/ez2
461 462 #------------------------------------------------------------------------------
463 -def OmegaRz(z):
464 """ Calculates the radiation density of the universe at redshift z. 465 466 @type z: float 467 @param z: redshift 468 @rtype: float 469 @return: radiation density of universe at redshift z 470 471 """ 472 ez2 = Ez2(z) 473 474 return OMEGA_R0*math.pow(1+z, 4)/ez2
475 476 #------------------------------------------------------------------------------
477 -def DeltaVz(z):
478 """Calculates the density contrast of a virialised region S{Delta}V(z), 479 assuming a S{Lambda}CDM-type flat cosmology. See, e.g., Bryan & Norman 480 1998 (ApJ, 495, 80). 481 482 @type z: float 483 @param z: redshift 484 @rtype: float 485 @return: density contrast of a virialised region at redshift z 486 487 @note: If OMEGA_M0+OMEGA_L0 is not equal to 1, this routine exits and 488 prints an error 489 message to the console. 490 491 """ 492 493 OMEGA_K = 1.0 - OMEGA_M0 - OMEGA_L0 494 495 if OMEGA_K == 0.0: 496 Omega_Mz = OmegaMz(z) 497 deltaVz = (18.0*math.pow(math.pi, 2)+82.0*(Omega_Mz-1.0)-39.0 * 498 math.pow(Omega_Mz-1, 2)) 499 return deltaVz 500 else: 501 raise Exception("cosmology is NOT flat.")
502 503 #------------------------------------------------------------------------------
504 -def RVirialXRayCluster(kT, z, betaT):
505 """Calculates the virial radius (in Mpc) of a galaxy cluster at redshift z 506 with X-ray temperature kT, assuming self-similar evolution and a flat 507 cosmology. See Arnaud et al. 2002 (A&A, 389, 1) and Bryan & Norman 1998 508 (ApJ, 495, 80). A flat S{Lambda}CDM-type flat cosmology is assumed. 509 510 @type kT: float 511 @param kT: cluster X-ray temperature in keV 512 @type z: float 513 @param z: redshift 514 @type betaT: float 515 @param betaT: the normalisation of the virial relation, for which Evrard et 516 al. 1996 (ApJ,469, 494) find a value of 1.05 517 @rtype: float 518 @return: virial radius of cluster in Mpc 519 520 @note: If OMEGA_M0+OMEGA_L0 is not equal to 1, this routine exits and 521 prints an error message to the console. 522 523 """ 524 525 OMEGA_K = 1.0 - OMEGA_M0 - OMEGA_L0 526 527 if OMEGA_K == 0.0: 528 Omega_Mz = OmegaMz(z) 529 deltaVz = (18.0 * math.pow(math.pi, 2) + 82.0 * (Omega_Mz-1.0)- 39.0 * 530 math.pow(Omega_Mz-1, 2)) 531 deltaz = (deltaVz*OMEGA_M0)/(18.0*math.pow(math.pi, 2)*Omega_Mz) 532 533 # The equation quoted in Arnaud, Aghanim & Neumann is for h50, so need 534 # to scale it 535 h50 = H0/50.0 536 Rv = (3.80*math.sqrt(betaT)*math.pow(deltaz, -0.5) * 537 math.pow(1.0+z, (-3.0/2.0)) * math.sqrt(kT/10.0)*(1.0/h50)) 538 539 return Rv 540 541 else: 542 raise Exception("cosmology is NOT flat.")
543 544 #------------------------------------------------------------------------------ 545