#+ # VIDE -- Void IDentification and Examination -- ./python_tools/vide/apTools/chi2/cosmologyTools.py # Copyright (C) 2010-2014 Guilhem Lavaux # Copyright (C) 2011-2014 P. M. Sutter # # This program is free software; you can redistribute it and/or modify # it under the terms of the GNU General Public License as published by # the Free Software Foundation; version 2 of the License. # # # This program is distributed in the hope that it will be useful, # but WITHOUT ANY WARRANTY; without even the implied warranty of # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the # GNU General Public License for more details. # # You should have received a copy of the GNU General Public License along # with this program; if not, write to the Free Software Foundation, Inc., # 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. #+ # a suite of functions to compute expansion rates, angular diameter # distances, and expected void stretching import numpy as np import scipy.integrate as integrate import healpy as healpy import os from vide.backend import * __all__=['expansion', 'angularDiameter', 'aveExpansion', 'getSurveyProps'] # returns 1/E(z) for the given cosmology def expansion(z, Om = 0.27, Ot = 1.0, w0 = -1.0, wa = 0.0): ez = Om * (1+z)**3 + (Ot-Om)# * (1+z)**(3.+3*wz) #ez = Om * (1+z)**3 + (Ot-Om)# * integrade.quad(eosDE, 0.0, z, args=(w0,wa))[0] ez = 1./np.sqrt(ez) return ez # returns DE value at redshift z def eosDE(z, w0 = -1.0, wa = 0.0): return w0 + wa*z/(1+z) # returns D_A(z) for the given cosmology def angularDiameter(z, Om = 0.27, Ot = 1.0, w0 = -1.0, wa = 0.0): da = integrate.quad(expansion, 0.0, z, args=(Om, Ot, w0, wa))[0] return da # ----------------------------------------------------------------------------- # returns average expected expansion for a given redshift range def aveExpansion(zStart, zEnd, Om = 0.27, Ot = 1.0, w0 = -1.0, wa = 0.0): if zStart == 0.0: zStart = 1.e-6 ave = integrate.quad(expansion, zStart, zEnd, args=(Om, Ot, w0, wa))[0] ave = (zEnd-zStart)/ave return ave # ----------------------------------------------------------------------------- # returns the volume and galaxy density for a given redshit slice def getSurveyProps(maskFile, zmin, zmax, selFunMin, selFunMax, portion, selectionFuncFile=None, useComoving=False): LIGHT_SPEED = 299792.458 mask = healpy.read_map(maskFile) area = (1.*np.size(np.where(mask > 0)) / np.size(mask)) * 4.*np.pi if useComoving: zmin = LIGHT_SPEED/100.*angularDiameter(zmin, Om=0.27) zmax = LIGHT_SPEED/100.*angularDiameter(zmax, Om=0.27) selFunMin = LIGHT_SPEED/100.*angularDiameter(selFunMin, Om=0.27) selFunMax = LIGHT_SPEED/100.*angularDiameter(selFunMax, Om=0.27) else: zmin = zmin * 3000 zmax = zmax * 3000 selFunMin *= 3000 selFunMax *= 3000 volume = area * (zmax**3 - zmin**3) / 3 if selectionFuncFile == None: nbar = 1.0 elif not os.access(selectionFuncFile, os.F_OK): print(" Warning, selection function file %s not found, using default of uniform selection." % selectionFuncFile) nbar = 1.0 else: selfunc = np.genfromtxt(selectionFuncFile) selfunc = np.array(selfunc) selfunc[:,0] = selfunc[:,0]/100. selfuncUnity = selfunc selfuncUnity[:,1] = 1.0 selfuncMin = selfunc[0,0] selfuncMax = selfunc[-1,0] selfuncDx = selfunc[1,0] - selfunc[0,0] selfuncN = np.size(selfunc[:,0]) selFunMin = max(selFunMin, selfuncMin) selFunMax = min(selFunMax, selfuncMax) def f(z): return selfunc[np.ceil((z-selfuncMin)/selfuncDx), 1]*z**2 def fTotal(z): return selfuncUnity[np.ceil((z-selfuncMin)/selfuncDx), 1]*z**2 zrange = np.linspace(selFunMin, selFunMax) nbar = scipy.integrate.quad(f, selFunMin, selFunMax) nbar = nbar[0] ntotal = scipy.integrate.quad(fTotal, 0.0, max(selfuncUnity[:,0])) ntotal = ntotal[0] nbar = ntotal / area / nbar return (volume, nbar)