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

Source Code for Module astLib.astCalc

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