1
2 """module for performing common calculations
3
4 (c) 2007-2011 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
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
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
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
81 N=1000
82
83
84 xMax=1.0
85 xMin=1.0/(1.0+z)
86 xStep=float((xMax-xMin)/N)
87
88
89 yTotal=0
90
91
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
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:
103 yTotal=yTotal+(4*yn)
104 else:
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
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
134 N=1000
135
136
137 xMax=1.0
138 xMin=1.0/(1.0+z)
139 xStep=float((xMax-xMin)/N)
140
141
142 yTotal=0
143
144
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
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:
156 yTotal=yTotal+(4*yn)
157 else:
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
168 """Calculates the line of sight comoving volume element per steradian dV/dz at redshift z.
169
170 @type z: float
171 @param z: redshift
172 @rtype: float
173 @return: comoving volume element per steradian
174
175 """
176
177 dH=C_LIGHT/H0
178 dVcdz=(dH*(da(z)**2)*((1+z)**2))/Ez(z)
179
180 return dVcdz
181
182
183 -def dl2z(distanceMpc):
184 """Calculates the redshift z corresponding to the luminosity distance given in Mpc.
185
186 @type distanceMpc: float
187 @param distanceMpc: distance in Mpc
188 @rtype: float
189 @return: redshift
190
191 """
192
193 dTarget=distanceMpc
194
195 toleranceMpc=0.1
196
197 zMin=0.0
198 zMax=10.0
199
200 diff=dl(zMax)-dTarget
201 while diff<0:
202 zMax=zMax+5.0
203 diff=dl(zMax)-dTarget
204
205 zTrial=zMin+(zMax-zMin)/2.0
206
207 dTrial=dl(zTrial)
208 diff=dTrial-dTarget
209 while abs(diff)>toleranceMpc:
210
211 if diff>0:
212 zMax=zMax-(zMax-zMin)/2.0
213 else:
214 zMin=zMin+(zMax-zMin)/2.0
215
216 zTrial=zMin+(zMax-zMin)/2.0
217 dTrial=dl(zTrial)
218 diff=dTrial-dTarget
219
220 return zTrial
221
222
223 -def dc2z(distanceMpc):
224 """Calculates the redshift z corresponding to the comoving distance given in Mpc.
225
226 @type distanceMpc: float
227 @param distanceMpc: distance in Mpc
228 @rtype: float
229 @return: redshift
230
231 """
232
233 dTarget=distanceMpc
234
235 toleranceMpc=0.1
236
237 zMin=0.0
238 zMax=10.0
239
240 diff=dc(zMax)-dTarget
241 while diff<0:
242 zMax=zMax+5.0
243 diff=dc(zMax)-dTarget
244
245 zTrial=zMin+(zMax-zMin)/2.0
246
247 dTrial=dc(zTrial)
248 diff=dTrial-dTarget
249 while abs(diff)>toleranceMpc:
250
251 if diff>0:
252 zMax=zMax-(zMax-zMin)/2.0
253 else:
254 zMin=zMin+(zMax-zMin)/2.0
255
256 zTrial=zMin+(zMax-zMin)/2.0
257 dTrial=dc(zTrial)
258 diff=dTrial-dTarget
259
260 return zTrial
261
262
264 """Calculates the age of the universe in Gyr at z=0 for the current set of cosmological
265 parameters.
266
267 @rtype: float
268 @return: age of the universe in Gyr at z=0
269
270 """
271
272 OMEGA_K=1.0-OMEGA_M0-OMEGA_L
273
274
275 N=1000
276
277
278 xMax=1.0
279 xMin=0.0
280 xStep=float((xMax-xMin)/N)
281
282
283 yTotal=0
284
285
286 for n in range(1,N):
287 x=xMin+xStep*n
288 yn=x/math.sqrt(OMEGA_M0*x+OMEGA_L*math.pow(x, 4)+OMEGA_K*math.pow(x, 2))
289
290
291 if n==0 or n==N:
292 yTotal=yTotal+yn
293 else:
294 oddness=float(n/2.0);
295 fractPart=math.modf(oddness)[0]
296 if fractPart==0.5:
297 yTotal=yTotal+(4*yn)
298 else:
299 yTotal=yTotal+(2*yn)
300
301 integralValue=(xMax-xMin)/(3*N)*yTotal
302
303 T0=(1.0/H0*integralValue*3.08e19)/3.16e7/1e9
304
305 return T0
306
307
309 """ Calculates the lookback time in Gyr to redshift z for the current set of cosmological
310 parameters.
311
312 @type z: float
313 @param z: redshift
314 @rtype: float
315 @return: lookback time in Gyr to redshift z
316
317 """
318 OMEGA_K=1.0-OMEGA_M0-OMEGA_L
319
320
321 N=1000
322
323
324 xMax=1.0
325 xMin=1.0/(1.0+z)
326 xStep=float((xMax-xMin)/N)
327
328
329 yTotal=0
330
331
332 for n in range(1,N):
333 x=xMin+xStep*n
334 yn=x/math.sqrt(OMEGA_M0*x+OMEGA_L*math.pow(x, 4)+OMEGA_K*math.pow(x, 2))
335
336
337 if n==0 or n==N:
338 yTotal=yTotal+yn
339 else:
340 oddness=float(n/2.0);
341 fractPart=math.modf(oddness)[0]
342 if fractPart==0.5:
343 yTotal=yTotal+(4*yn)
344 else:
345 yTotal=yTotal+(2*yn)
346
347 integralValue=(xMax-xMin)/(3*N)*yTotal
348
349 T0=(1.0/H0*integralValue*3.08e19)/3.16e7/1e9
350
351 return T0
352
353
355 """Calculates the age of the universe at redshift z for the current set of cosmological
356 parameters.
357
358 @type z: float
359 @param z: redshift
360 @rtype: float
361 @return: age of the universe in Gyr at redshift z
362
363 """
364
365 TZ=t0()-tl(z)
366
367 return TZ
368
369
371 """Calculates the redshift z corresponding to lookback time tlGyr given in Gyr.
372
373 @type tlGyr: float
374 @param tlGyr: lookback time in Gyr
375 @rtype: float
376 @return: redshift
377
378 """
379
380 tTarget=tlGyr
381
382 toleranceGyr=0.001
383
384 zMin=0.0
385 zMax=10.0
386
387 diff=tl(zMax)-tTarget
388 while diff<0:
389 zMax=zMax+5.0
390 diff=tl(zMax)-tTarget
391
392 zTrial=zMin+(zMax-zMin)/2.0
393
394 tTrial=tl(zTrial)
395 diff=tTrial-tTarget
396 while abs(diff)>toleranceGyr:
397
398 if diff>0:
399 zMax=zMax-(zMax-zMin)/2.0
400 else:
401 zMin=zMin+(zMax-zMin)/2.0
402
403 zTrial=zMin+(zMax-zMin)/2.0
404 tTrial=tl(zTrial)
405 diff=tTrial-tTarget
406
407 return zTrial
408
409
411 """Calculates the redshift z corresponding to age of the universe tzGyr given in Gyr.
412
413 @type tzGyr: float
414 @param tzGyr: age of the universe in Gyr
415 @rtype: float
416 @return: redshift
417
418 """
419 tl=t0()-tzGyr
420 z=tl2z(tl)
421
422 return z
423
424
426 """Calculates the absolute magnitude of an object at given luminosity distance in Mpc.
427
428 @type appMag: float
429 @param appMag: apparent magnitude of object
430 @type distMpc: float
431 @param distMpc: distance to object in Mpc
432 @rtype: float
433 @return: absolute magnitude of object
434
435 """
436 absMag=appMag-(5.0*math.log10(distMpc*1.0e5))
437
438 return absMag
439
440
442 """Calculates the value of E(z), which describes evolution of the Hubble parameter with
443 redshift, at redshift z for the current set of cosmological parameters. See, e.g., Bryan &
444 Norman 1998 (ApJ, 495, 80).
445
446 @type z: float
447 @param z: redshift
448 @rtype: float
449 @return: value of E(z) at redshift z
450
451 """
452
453 Ez=math.sqrt(Ez2(z))
454
455 return Ez
456
457
459 """Calculates the value of E(z)^2, which describes evolution of the Hubble parameter with
460 redshift, at redshift z for the current set of cosmological parameters. See, e.g., Bryan &
461 Norman 1998 (ApJ, 495, 80).
462
463 @type z: float
464 @param z: redshift
465 @rtype: float
466 @return: value of E(z)^2 at redshift z
467
468 """
469
470 Ez2=OMEGA_M0*math.pow(1.0+z, 3)+(1.0-OMEGA_M0-OMEGA_L)*math.pow(1.0+z, 2)+OMEGA_L
471
472 return Ez2
473
474
476 """Calculates the matter density of the universe at redshift z. See, e.g., Bryan & Norman
477 1998 (ApJ, 495, 80).
478
479 @type z: float
480 @param z: redshift
481 @rtype: float
482 @return: matter density of universe at redshift z
483
484 """
485 ez2=Ez2(z)
486
487 Omega_Mz=(OMEGA_M0*math.pow(1.0+z, 3))/ez2
488
489 return Omega_Mz
490
491
493 """Calculates the density contrast of a virialised region S{Delta}V(z), assuming a
494 S{Lambda}CDM-type flat cosmology. See, e.g., Bryan & Norman 1998 (ApJ, 495, 80).
495
496 @type z: float
497 @param z: redshift
498 @rtype: float
499 @return: density contrast of a virialised region at redshift z
500
501 @note: If OMEGA_M0+OMEGA_L is not equal to 1, this routine exits and prints an error
502 message to the console.
503
504 """
505
506 OMEGA_K=1.0-OMEGA_M0-OMEGA_L
507
508 if OMEGA_K==0.0:
509 Omega_Mz=OmegaMz(z)
510 deltaVz=18.0*math.pow(math.pi, 2)+82.0*(Omega_Mz-1.0)-39.0*math.pow(Omega_Mz-1, 2)
511 return deltaVz
512 else:
513 raise Exception, "cosmology is NOT flat."
514
515
517 """Calculates the virial radius (in Mpc) of a galaxy cluster at redshift z with X-ray temperature
518 kT, assuming self-similar evolution and a flat cosmology. See Arnaud et al. 2002 (A&A,
519 389, 1) and Bryan & Norman 1998 (ApJ, 495, 80). A flat S{Lambda}CDM-type flat cosmology is
520 assumed.
521
522 @type kT: float
523 @param kT: cluster X-ray temperature in keV
524 @type z: float
525 @param z: redshift
526 @type betaT: float
527 @param betaT: the normalisation of the virial relation, for which Evrard et al. 1996
528 (ApJ,469, 494) find a value of 1.05
529 @rtype: float
530 @return: virial radius of cluster in Mpc
531
532 @note: If OMEGA_M0+OMEGA_L is not equal to 1, this routine exits and prints an error
533 message to the console.
534
535 """
536
537 OMEGA_K=1.0-OMEGA_M0-OMEGA_L
538
539 if OMEGA_K==0.0:
540 Omega_Mz=OmegaMz(z)
541 deltaVz=18.0*math.pow(math.pi, 2)+82.0*(Omega_Mz-1.0)-39.0*math.pow(Omega_Mz-1, 2)
542 deltaz=(deltaVz*OMEGA_M0)/(18.0*math.pow(math.pi, 2)*Omega_Mz)
543
544
545 h50=H0/50.0
546 Rv=3.80*math.sqrt(betaT)*math.pow(deltaz, -0.5) \
547 *math.pow(1.0+z, (-3.0/2.0))*math.sqrt(kT/10.0)*(1.0/h50)
548
549 return Rv
550
551 else:
552 raise Exception, "cosmology is NOT flat."
553
554
555